CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/src/Calibration/EcalCalibAlgos/src/PhiSymmetryCalibration.cc

Go to the documentation of this file.
00001 #include "Calibration/EcalCalibAlgos/interface/PhiSymmetryCalibration.h"
00002 
00003 // System include files
00004 #include <memory>
00005 
00006 // Framework
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 #include "FWCore/Framework/interface/ESHandle.h"
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 
00011 #include "DataFormats/EcalRecHit/interface/EcalRecHitCollections.h"
00012 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00013 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00014 #include "CondFormats/DataRecord/interface/EcalIntercalibConstantsRcd.h"
00015 #include "CondFormats/EcalObjects/interface/EcalIntercalibErrors.h"
00016 #include "CondTools/Ecal/interface/EcalIntercalibConstantsXMLTranslator.h"
00017 #include "FWCore/Framework/interface/LuminosityBlock.h"
00018 
00019 // Geometry
00020 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00021 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00022 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00023 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00024 
00025 //Channel status
00026 #include "CondFormats/EcalObjects/interface/EcalChannelStatus.h"
00027 #include "CondFormats/DataRecord/interface/EcalChannelStatusRcd.h"
00028 #include "CondFormats/EcalObjects/interface/EcalChannelStatusCode.h"
00029 
00030 #include "FWCore/Framework/interface/Run.h"
00031 
00032 
00033 
00034 
00035 
00036 
00037 #include "boost/filesystem/operations.hpp"
00038 
00039 using namespace std;
00040 #include <fstream>
00041 #include <iostream>
00042 #include "TH2F.h"
00043 #include "TFile.h"
00044 #include "TTree.h"
00045 #include "TH1F.h"
00046 #include "TF1.h"
00047 #include "TGraph.h"
00048 #include "TCanvas.h"
00049 
00050 const float PhiSymmetryCalibration::kMiscalRangeEB = .05;
00051 const float PhiSymmetryCalibration::kMiscalRangeEE = .10;
00052 
00053 
00054 
00055 
00056 //_____________________________________________________________________________
00057 // Class constructor
00058 
00059 PhiSymmetryCalibration::PhiSymmetryCalibration(const edm::ParameterSet& iConfig) :
00060 
00061   ecalHitsProducer_(iConfig.getParameter<std::string>("ecalRecHitsProducer")),
00062   barrelHits_( iConfig.getParameter< std::string > ("barrelHitCollection")),
00063   endcapHits_( iConfig.getParameter< std::string > ("endcapHitCollection")),
00064   eCut_barl_( iConfig.getParameter< double > ("eCut_barrel") ),
00065   ap_( iConfig.getParameter<double> ("ap") ),
00066   b_( iConfig.getParameter<double> ("b") ), 
00067   eventSet_( iConfig.getParameter< int > ("eventSet") ),
00068   statusThreshold_(iConfig.getUntrackedParameter<int>("statusThreshold",3)),
00069   reiteration_(iConfig.getUntrackedParameter< bool > ("reiteration",false)),
00070   oldcalibfile_(iConfig.getUntrackedParameter<std::string>("oldcalibfile",
00071                                             "EcalintercalibConstants.xml"))
00072 {
00073 
00074 
00075   isfirstpass_=true;
00076 
00077   et_spectrum_b_histos.resize(kBarlRings);
00078   e_spectrum_b_histos.resize(kBarlRings);
00079   et_spectrum_e_histos.resize(kEndcEtaRings);
00080   e_spectrum_e_histos.resize(kEndcEtaRings); 
00081   
00082   spectra=true;
00083 
00084   nevents_=0;
00085   eventsinrun_=0;
00086   eventsinlb_=0;
00087 }
00088 
00089 
00090 //_____________________________________________________________________________
00091 // Close files, etc.
00092 
00093 PhiSymmetryCalibration::~PhiSymmetryCalibration()
00094 {
00095 
00096 
00097   for(Int_t i=0;i<kBarlRings;i++){
00098     delete et_spectrum_b_histos[i];
00099     delete e_spectrum_b_histos[i];
00100     
00101   }
00102   for(Int_t i=0;i<kEndcEtaRings;i++){
00103     delete et_spectrum_e_histos[i];
00104     delete e_spectrum_e_histos[i];    
00105   }
00106 
00107   
00108 }
00109 
00110 
00111 //_____________________________________________________________________________
00112 // Initialize algorithm
00113 
00114 void PhiSymmetryCalibration::beginJob( )
00115 {
00116 
00117 
00118   // initialize arrays
00119   for (int sign=0; sign<kSides; sign++) {
00120     for (int ieta=0; ieta<kBarlRings; ieta++) {
00121       for (int iphi=0; iphi<kBarlWedges; iphi++) {
00122         etsum_barl_[ieta][iphi][sign]=0.;
00123         nhits_barl_[ieta][iphi][sign]=0;
00124 
00125       }
00126     }
00127     for (int ix=0; ix<kEndcWedgesX; ix++) {
00128       for (int iy=0; iy<kEndcWedgesY; iy++) {
00129         etsum_endc_[ix][iy][sign]=0.;
00130         nhits_endc_[ix][iy][sign]=0;
00131       }
00132     }
00133   }
00134 
00135 
00136 
00137   for (int imiscal=0; imiscal<kNMiscalBinsEB; imiscal++) {
00138     miscalEB_[imiscal]= (1-kMiscalRangeEB) + float(imiscal)* (2*kMiscalRangeEB/(kNMiscalBinsEB-1));
00139     for (int ieta=0; ieta<kBarlRings; ieta++) etsum_barl_miscal_[imiscal][ieta]=0.;
00140   }
00141 
00142   for (int imiscal=0; imiscal<kNMiscalBinsEE; imiscal++) {
00143     miscalEE_[imiscal]= (1-kMiscalRangeEE) + float(imiscal)* (2*kMiscalRangeEE/(kNMiscalBinsEE-1));
00144     for (int ring=0; ring<kEndcEtaRings; ring++) etsum_endc_miscal_[imiscal][ring]=0.;
00145   }
00146 
00147 
00148 
00149   // start spectra stuff
00150   if (eventSet_!=1) spectra = false;
00151   
00152   if(spectra)
00153     {
00154       ostringstream t;
00155       for(Int_t i=0;i<kBarlRings;i++)
00156         {
00157           t << "et_spectrum_b_" << i+1;
00158           et_spectrum_b_histos[i]=new TH1F(t.str().c_str(),";E_{T} [MeV]",50,0.,500.);
00159           t.str("");
00160           
00161           t << "e_spectrum_b_" << i+1;
00162           e_spectrum_b_histos[i]=new TH1F(t.str().c_str(),";E [MeV]",50,0.,500.);
00163           t.str("");
00164         
00165         }
00166       for(Int_t i=0;i<kEndcEtaRings;i++)
00167         {
00168           t << "et_spectrum_e_" << i+1;
00169           et_spectrum_e_histos[i]=new TH1F(t.str().c_str(),";E_{T} [MeV]",75,0.,1500.);
00170           t.str("");
00171           
00172           t << "e_spectrum_e_" << i+1;
00173           e_spectrum_e_histos[i]=new TH1F(t.str().c_str(),";E [MeV]",75,0.,1500.);
00174           t.str("");
00175         
00176         }
00177     }
00178   // end spectra stuff
00179 }
00180 
00181 
00182 //_____________________________________________________________________________
00183 // Terminate algorithm
00184 
00185 void PhiSymmetryCalibration::endJob()
00186 {
00187 
00188   edm::LogInfo("Calibration") << "[PhiSymmetryCalibration] At end of job";
00189 
00190   // start spectra stuff
00191   if(spectra)
00192     {
00193       TFile f("Espectra_plus.root","recreate");
00194 
00195       for(int i=0;i<kBarlRings;i++){
00196         et_spectrum_b_histos[i]->Write();
00197         e_spectrum_b_histos[i]->Write();
00198       }
00199 
00200       for(int i=0;i<kEndcEtaRings;i++){
00201         et_spectrum_e_histos[i]->Write();
00202         e_spectrum_e_histos[i]->Write();
00203       }
00204 
00205       f.Close();
00206     }
00207   
00208 
00209 
00210 
00211   if (eventSet_==1) {
00212     // calculate factors to convert from fractional deviation of ET sum from 
00213     // the mean to the estimate of the miscalibration factor
00214     getKfactors();
00215 
00216     std::ofstream k_barl_out("k_barl.dat", ios::out);
00217     for (int ieta=0; ieta<kBarlRings; ieta++)
00218       k_barl_out << ieta << " " << k_barl_[ieta] << endl;
00219     k_barl_out.close();
00220 
00221     std::ofstream k_endc_out("k_endc.dat", ios::out);
00222     for (int ring=0; ring<kEndcEtaRings; ring++)
00223       k_endc_out << ring << " " << k_endc_[ring] << endl;
00224     k_endc_out.close();
00225   }
00226 
00227 
00228   if (eventSet_!=0) {
00229     // output ET sums
00230 
00231     stringstream etsum_file_barl;
00232     etsum_file_barl << "etsum_barl_"<<eventSet_<<".dat";
00233 
00234     std::ofstream etsum_barl_out(etsum_file_barl.str().c_str(),ios::out);
00235 
00236     for (int ieta=0; ieta<kBarlRings; ieta++) {
00237       for (int iphi=0; iphi<kBarlWedges; iphi++) {
00238         for (int sign=0; sign<kSides; sign++) {
00239           etsum_barl_out << eventSet_ << " " << ieta << " " << iphi << " " << sign 
00240                          << " " << etsum_barl_[ieta][iphi][sign] << " " 
00241                          <<  nhits_barl_[ieta][iphi][sign] << endl;
00242         }
00243       }
00244     }
00245     etsum_barl_out.close();
00246 
00247     stringstream etsum_file_endc;
00248     etsum_file_endc << "etsum_endc_"<<eventSet_<<".dat";
00249 
00250     std::ofstream etsum_endc_out(etsum_file_endc.str().c_str(),ios::out);
00251     for (int ix=0; ix<kEndcWedgesX; ix++) {
00252       for (int iy=0; iy<kEndcWedgesY; iy++) {
00253         int ring = e_.endcapRing_[ix][iy];
00254         if (ring!=-1) {
00255           for (int sign=0; sign<kSides; sign++) {
00256             etsum_endc_out << eventSet_ << " " << ix << " " << iy << " " << sign 
00257                            << " " << etsum_endc_[ix][iy][sign] << " " 
00258                            << nhits_endc_[ix][iy][sign]<<" " 
00259                            << e_.endcapRing_[ix][iy]<<endl;
00260           }
00261         }
00262       }
00263     }
00264     etsum_endc_out.close();
00265   } 
00266   cout<<"Events processed " << nevents_<< endl;
00267 }
00268 
00269 
00270 //_____________________________________________________________________________
00271 // Called at each event
00272 
00273 void PhiSymmetryCalibration::analyze( const edm::Event& event, const edm::EventSetup& setup )
00274 {
00275   using namespace edm;
00276   using namespace std;
00277 
00278   if (isfirstpass_) {
00279     setUp(setup);
00280     isfirstpass_=false;
00281   }
00282 
00283   
00284   Handle<EBRecHitCollection> barrelRecHitsHandle;
00285   Handle<EERecHitCollection> endcapRecHitsHandle;
00286   
00287   event.getByLabel(ecalHitsProducer_,barrelHits_,barrelRecHitsHandle);
00288   if (!barrelRecHitsHandle.isValid()) {
00289     LogError("") << "[PhiSymmetryCalibration] Error! Can't get product!" << std::endl;
00290   }
00291   
00292   event.getByLabel(ecalHitsProducer_,endcapHits_,endcapRecHitsHandle);
00293   if (!endcapRecHitsHandle.isValid()) {
00294     LogError("") << "[PhiSymmetryCalibration] Error! Can't get product!" << std::endl;
00295   }
00296   
00297  
00298   // get the ecal geometry
00299   edm::ESHandle<CaloGeometry> geoHandle;
00300   setup.get<CaloGeometryRecord>().get(geoHandle);
00301   const CaloSubdetectorGeometry *barrelGeometry = 
00302     geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalBarrel);
00303   const CaloSubdetectorGeometry *endcapGeometry = 
00304     geoHandle->getSubdetectorGeometry(DetId::Ecal, EcalEndcap);
00305  
00306   bool pass=false;
00307   // select interesting EcalRecHits (barrel)
00308   EBRecHitCollection::const_iterator itb;
00309   for (itb=barrelRecHitsHandle->begin(); itb!=barrelRecHitsHandle->end(); itb++) {
00310     EBDetId hit = EBDetId(itb->id());
00311     float eta = barrelGeometry->getGeometry(hit)->getPosition().eta();
00312     float et = itb->energy()/cosh(eta);
00313     float e  = itb->energy();
00314     
00315     
00316 
00317     // if iterating, correct by the previous calib constants found,
00318     // which are supplied in the form of correction 
00319     if (reiteration_) {
00320       et= et  * oldCalibs_[hit];
00321       e = e  * oldCalibs_[hit];
00322     }
00323 
00324     float et_thr = eCut_barl_/cosh(eta) + 1.;
00325 
00326     int sign = hit.ieta()>0 ? 1 : 0;
00327 
00328     if (e >  eCut_barl_ && et < et_thr && e_.goodCell_barl[abs(hit.ieta())-1][hit.iphi()-1][sign]) {
00329       etsum_barl_[abs(hit.ieta())-1][hit.iphi()-1][sign] += et;
00330       nhits_barl_[abs(hit.ieta())-1][hit.iphi()-1][sign] ++;
00331       pass =true;
00332     }//if energy
00333 
00334     if (eventSet_==1) {
00335       // apply a miscalibration to all crystals and increment the 
00336       // ET sum, combined for all crystals
00337       for (int imiscal=0; imiscal<kNMiscalBinsEB; imiscal++) {
00338         if (miscalEB_[imiscal]*e >  eCut_barl_&& miscalEB_[imiscal]*et < et_thr && e_.goodCell_barl[abs(hit.ieta())-1][hit.iphi()-1][sign]) {
00339           etsum_barl_miscal_[imiscal][abs(hit.ieta())-1] += miscalEB_[imiscal]*et;
00340         }
00341       }
00342 
00343       // spectra stuff
00344       if(spectra && hit.ieta()>0) //POSITIVE!!!
00345         //      if(spectra && hit.ieta()<0) //NEGATIVE!!!
00346         {
00347           et_spectrum_b_histos[abs(hit.ieta())-1]->Fill(et*1000.);
00348           e_spectrum_b_histos[abs(hit.ieta())-1]->Fill(e*1000.);
00349         }//if spectra
00350       
00351     }//if eventSet_==1
00352   }//for barl
00353 
00354 
00355   // select interesting EcalRecHits (endcaps)
00356   EERecHitCollection::const_iterator ite;
00357   for (ite=endcapRecHitsHandle->begin(); ite!=endcapRecHitsHandle->end(); ite++) {
00358     EEDetId hit = EEDetId(ite->id());
00359     float eta = abs(endcapGeometry->getGeometry(hit)->getPosition().eta());
00360     //float phi = endcapGeometry->getGeometry(hit)->getPosition().phi();
00361 
00362     float et = ite->energy()/cosh(eta);
00363     float e  = ite->energy();
00364 
00365     // if iterating, multiply by the previous correction factor
00366     if (reiteration_) {
00367       et= et * oldCalibs_[hit];
00368       e = e * oldCalibs_[hit];
00369     }
00370 
00371     int sign = hit.zside()>0 ? 1 : 0;
00372 
00373 
00374     // changes of eCut_endc_ -> variable linearthr 
00375     // e_cut = ap + eta_ring*b
00376 
00377     double eCut_endc=0;
00378     for (int ring=0; ring<kEndcEtaRings; ring++) {
00379 
00380       if(eta>e_.etaBoundary_[ring] && eta<e_.etaBoundary_[ring+1])
00381         {  
00382           float eta_ring= abs(e_.cellPos_[ring][50].eta())  ;
00383           eCut_endc = ap_ + eta_ring*b_;
00384 
00385         }
00386     }
00387 
00388 
00389     float et_thr = eCut_endc/cosh(eta) + 1.;
00390    
00391     if (e > eCut_endc && et < et_thr && e_.goodCell_endc[hit.ix()-1][hit.iy()-1][sign]){
00392       etsum_endc_[hit.ix()-1][hit.iy()-1][sign] += et;
00393       nhits_endc_[hit.ix()-1][hit.iy()-1][sign] ++;
00394       pass=true;
00395     }
00396  
00397    
00398 
00399     if (eventSet_==1) {
00400       // apply a miscalibration to all crystals and increment the 
00401       // ET sum, combined for all crystals
00402       for (int imiscal=0; imiscal<kNMiscalBinsEE; imiscal++) {
00403         if (miscalEE_[imiscal]*e> eCut_endc && et*miscalEE_[imiscal] < et_thr && e_.goodCell_endc[hit.ix()-1][hit.iy()-1][sign]){
00404           int ring = e_.endcapRing_[hit.ix()-1][hit.iy()-1];
00405           etsum_endc_miscal_[imiscal][ring] += miscalEE_[imiscal]*et;
00406         }
00407       }
00408 
00409       // spectra stuff
00410       if(spectra && hit.zside()>0) //POSITIVE!!!
00411 
00412         {
00413           int ring = e_.endcapRing_[hit.ix()-1][hit.iy()-1];
00414 
00415           et_spectrum_e_histos[ring]->Fill(et*1000.);
00416           e_spectrum_e_histos[ring]->Fill(e*1000.);
00417 
00418           if(ring==16)
00419             {
00420               //int iphi_endc = 0;
00421               for (int ip=0; ip<e_.nRing_[ring]; ip++) {
00422                 //if (phi==e_.phi_endc_[ip][ring]) iphi_endc=ip;
00423               }
00424 
00425             }
00426         }//if spectra
00427 
00428     }//if eventSet_==1
00429   }//for endc
00430 
00431   if (pass) {
00432     nevents_++;
00433     eventsinrun_++;
00434     eventsinlb_++;
00435   }
00436 }
00437 
00438 void PhiSymmetryCalibration::endRun(edm::Run& run, const edm::EventSetup&){
00439  
00440   
00441   std::cout  << "PHIREPRT : run "<< run.run() 
00442              << " start " << (run.beginTime().value()>>32)
00443              << " end "   << (run.endTime().value()>>32) 
00444              << " dur "   << (run.endTime().value()>>32)- (run.beginTime().value()>>32)
00445           
00446              << " npass "      << eventsinrun_  << std::endl;
00447   eventsinrun_=0;        
00448  
00449   return ;
00450 
00451 }
00452 
00453 //_____________________________________________________________________________
00454 
00455 void PhiSymmetryCalibration::getKfactors()
00456 {
00457 
00458   float epsilon_T_eb[kNMiscalBinsEB];
00459   float epsilon_M_eb[kNMiscalBinsEB];
00460 
00461   float epsilon_T_ee[kNMiscalBinsEE];
00462   float epsilon_M_ee[kNMiscalBinsEE];
00463 
00464   std::vector<TGraph*>  k_barl_graph(kBarlRings);
00465   std::vector<TCanvas*> k_barl_plot(kBarlRings);
00466 
00467   for (int ieta=0; ieta<kBarlRings; ieta++) {
00468     for (int imiscal=0; imiscal<kNMiscalBinsEB; imiscal++) {
00469       int middlebin =  int (kNMiscalBinsEB/2);
00470       epsilon_T_eb[imiscal] = etsum_barl_miscal_[imiscal][ieta]/etsum_barl_miscal_[middlebin][ieta] - 1.;
00471       epsilon_M_eb[imiscal] = miscalEB_[imiscal] - 1.;
00472     }
00473     k_barl_graph[ieta] = new TGraph (kNMiscalBinsEB,epsilon_M_eb,epsilon_T_eb);
00474     k_barl_graph[ieta]->Fit("pol1");
00475 
00476 
00477     ostringstream t;
00478     t<< "k_barl_" << ieta+1; 
00479     k_barl_plot[ieta] = new TCanvas(t.str().c_str(),"");
00480     k_barl_plot[ieta]->SetFillColor(10);
00481     k_barl_plot[ieta]->SetGrid();
00482     k_barl_graph[ieta]->SetMarkerSize(1.);
00483     k_barl_graph[ieta]->SetMarkerColor(4);
00484     k_barl_graph[ieta]->SetMarkerStyle(20);
00485     k_barl_graph[ieta]->GetXaxis()->SetLimits(-1.*kMiscalRangeEB,kMiscalRangeEB);
00486     k_barl_graph[ieta]->GetXaxis()->SetTitleSize(.05);
00487     k_barl_graph[ieta]->GetYaxis()->SetTitleSize(.05);
00488     k_barl_graph[ieta]->GetXaxis()->SetTitle("#epsilon_{M}");
00489     k_barl_graph[ieta]->GetYaxis()->SetTitle("#epsilon_{T}");
00490     k_barl_graph[ieta]->Draw("AP");
00491 
00492     k_barl_[ieta] = k_barl_graph[ieta]->GetFunction("pol1")->GetParameter(1);
00493     std::cout << "k_barl_[" << ieta << "]=" << k_barl_[ieta] << std::endl;
00494   }
00495 
00496 
00497   std::vector<TGraph*>  k_endc_graph(kEndcEtaRings);
00498   std::vector<TCanvas*> k_endc_plot(kEndcEtaRings);
00499 
00500   for (int ring=0; ring<kEndcEtaRings; ring++) {
00501     for (int imiscal=0; imiscal<kNMiscalBinsEE; imiscal++) {
00502       int middlebin =  int (kNMiscalBinsEE/2);
00503       epsilon_T_ee[imiscal] = etsum_endc_miscal_[imiscal][ring]/etsum_endc_miscal_[middlebin][ring] - 1.;
00504       epsilon_M_ee[imiscal] = miscalEE_[imiscal] - 1.;
00505     }
00506     k_endc_graph[ring] = new TGraph (kNMiscalBinsEE,epsilon_M_ee,epsilon_T_ee);
00507     k_endc_graph[ring]->Fit("pol1");
00508 
00509     ostringstream t;
00510     t<< "k_endc_"<< ring+1;
00511     k_endc_plot[ring] = new TCanvas(t.str().c_str(),"");
00512     k_endc_plot[ring]->SetFillColor(10);
00513     k_endc_plot[ring]->SetGrid();
00514     k_endc_graph[ring]->SetMarkerSize(1.);
00515     k_endc_graph[ring]->SetMarkerColor(4);
00516     k_endc_graph[ring]->SetMarkerStyle(20);
00517     k_endc_graph[ring]->GetXaxis()->SetLimits(-1*kMiscalRangeEE,kMiscalRangeEE);
00518     k_endc_graph[ring]->GetXaxis()->SetTitleSize(.05);
00519     k_endc_graph[ring]->GetYaxis()->SetTitleSize(.05);
00520     k_endc_graph[ring]->GetXaxis()->SetTitle("#epsilon_{M}");
00521     k_endc_graph[ring]->GetYaxis()->SetTitle("#epsilon_{T}");
00522     k_endc_graph[ring]->Draw("AP");
00523 
00524     k_endc_[ring] = k_endc_graph[ring]->GetFunction("pol1")->GetParameter(1);
00525     std::cout << "k_endc_[" << ring << "]=" << k_endc_[ring] << std::endl;
00526   }
00527  
00528   TFile f("PhiSymmetryCalibration_kFactors.root","recreate");
00529   for (int ieta=0; ieta<kBarlRings; ieta++) { 
00530     k_barl_plot[ieta]->Write();
00531     delete k_barl_plot[ieta]; 
00532     delete k_barl_graph[ieta];
00533   }
00534   for (int ring=0; ring<kEndcEtaRings; ring++) { 
00535     k_endc_plot[ring]->Write();
00536     delete k_endc_plot[ring];
00537     delete k_endc_graph[ring];
00538   }
00539   f.Close();
00540 
00541 }
00542 
00543 
00544 
00545 
00546 
00547 //_____________________________________________________________________________
00548 
00549 
00550 
00551 void PhiSymmetryCalibration::setUp(const edm::EventSetup& setup){
00552 
00553   edm::ESHandle<EcalChannelStatus> chStatus;
00554   setup.get<EcalChannelStatusRcd>().get(chStatus);
00555 
00556   edm::ESHandle<CaloGeometry> geoHandle;
00557   setup.get<CaloGeometryRecord>().get(geoHandle);
00558 
00559   e_.setup(&(*geoHandle), &(*chStatus), statusThreshold_);
00560  
00561   
00562   if (reiteration_){   
00563     
00564     EcalCondHeader h;
00565     // namespace fs = boost::filesystem;
00566 //     fs::path p(oldcalibfile_.c_str(),fs::native);
00567 //     if (!fs::exists(p)) edm::LogError("PhiSym") << "File not found: " 
00568 //                                              << oldcalibfile_ <<endl;
00569     
00570     edm::FileInPath fip("Calibration/EcalCalibAlgos/data/"+oldcalibfile_);
00571     
00572 
00573     
00574     int ret=
00575     EcalIntercalibConstantsXMLTranslator::readXML(fip.fullPath(),h,oldCalibs_);    
00576     if (ret) edm::LogError("PhiSym")<<"Error reading XML files"<<endl;;
00577     
00578   } else {
00579     // in fact if not reiterating, oldCalibs_ will never be used
00580     edm::ESHandle<EcalIntercalibConstants> pIcal;      
00581     setup.get<EcalIntercalibConstantsRcd>().get(pIcal);
00582     oldCalibs_=*pIcal;
00583 
00584   }
00585   
00586 }
00587 
00588 
00589 void PhiSymmetryCalibration::endLuminosityBlock(edm::LuminosityBlock const& lb, edm::EventSetup const&){
00590 
00591   
00592   if ((lb.endTime().value()>>32)- (lb.beginTime().value()>>32) <60 ) 
00593     return;
00594 
00595   std::cout  << "PHILB : run "<< lb.run()
00596              << " id " << lb.id() 
00597              << " start " << (lb.beginTime().value()>>32)
00598              << " end "   << (lb.endTime().value()>>32) 
00599              << " dur "   << (lb.endTime().value()>>32)- (lb.beginTime().value()>>32)
00600     
00601              << " npass "      << eventsinlb_  << std::endl;
00602   
00603   eventsinlb_=0;
00604 
00605 }