CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_3/src/RecoVertex/BeamSpotProducer/plugins/BeamSpotAnalyzer.cc

Go to the documentation of this file.
00001 
00015 // C++ standard
00016 #include <string>
00017 // CMS
00018 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00019 #include "RecoVertex/BeamSpotProducer/interface/BeamSpotAnalyzer.h"
00020 #include "RecoVertex/BeamSpotProducer/interface/BSFitter.h"
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include "FWCore/Framework/interface/LuminosityBlock.h"
00024 #include "FWCore/Framework/interface/MakerMacros.h"
00025 
00026 #include "TMath.h"
00027 
00028 BeamSpotAnalyzer::BeamSpotAnalyzer(const edm::ParameterSet& iConfig)
00029 {
00030   // get parameter
00031   write2DB_       = iConfig.getParameter<edm::ParameterSet>("BSAnalyzerParameters").getParameter<bool>("WriteToDB");
00032   runallfitters_  = iConfig.getParameter<edm::ParameterSet>("BSAnalyzerParameters").getParameter<bool>("RunAllFitters");
00033   fitNLumi_       = iConfig.getParameter<edm::ParameterSet>("BSAnalyzerParameters").getUntrackedParameter<int>("fitEveryNLumi",-1);
00034   resetFitNLumi_  = iConfig.getParameter<edm::ParameterSet>("BSAnalyzerParameters").getUntrackedParameter<int>("resetEveryNLumi",-1);
00035   runbeamwidthfit_  = iConfig.getParameter<edm::ParameterSet>("BSAnalyzerParameters").getParameter<bool>("RunBeamWidthFit");
00036 
00037   theBeamFitter = new BeamFitter(iConfig);
00038   theBeamFitter->resetTrkVector();
00039   theBeamFitter->resetLSRange();
00040   theBeamFitter->resetCutFlow();
00041   theBeamFitter->resetRefTime();
00042   theBeamFitter->resetPVFitter();
00043 
00044   ftotalevents = 0;
00045   ftmprun0 = ftmprun = -1;
00046   countLumi_ = 0;
00047   beginLumiOfBSFit_ = endLumiOfBSFit_ = -1;
00048   previousLumi_ = previousRun_ = 0;
00049   Org_resetFitNLumi_ = resetFitNLumi_;
00050 }
00051 
00052 
00053 BeamSpotAnalyzer::~BeamSpotAnalyzer()
00054 {
00055   delete theBeamFitter;
00056 }
00057 
00058 
00059 void
00060 BeamSpotAnalyzer::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup)
00061 {
00062         ftotalevents++;
00063         theBeamFitter->readEvent(iEvent);
00064         ftmprun = iEvent.id().run();
00065 
00066 }
00067 
00068 
00069 
00070 void
00071 BeamSpotAnalyzer::beginJob()
00072 {
00073 }
00074 
00075 //--------------------------------------------------------
00076 void
00077 BeamSpotAnalyzer::beginLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00078                                                                            const edm::EventSetup& context) {
00079 
00080     const edm::TimeValue_t fbegintimestamp = lumiSeg.beginTime().value();
00081     const std::time_t ftmptime = fbegintimestamp >> 32;
00082 
00083     if ( countLumi_ == 0 || (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) ) {
00084         ftmprun0 = lumiSeg.run();
00085         ftmprun = ftmprun0;
00086         beginLumiOfBSFit_ = lumiSeg.luminosityBlock();
00087         refBStime[0] = ftmptime;
00088     }
00089 
00090     countLumi_++;
00091     //std::cout << "Lumi # " << countLumi_ << std::endl;
00092     if ( ftmprun == previousRun_ ) {
00093       if ( (previousLumi_ + 1) != int(lumiSeg.luminosityBlock()) )
00094         edm::LogWarning("BeamSpotAnalyzer") << "LUMI SECTIONS ARE NOT SORTED!";
00095     }
00096 }
00097 
00098 //--------------------------------------------------------
00099 void
00100 BeamSpotAnalyzer::endLuminosityBlock(const edm::LuminosityBlock& lumiSeg,
00101                                                                          const edm::EventSetup& iSetup) {
00102 
00103   //LogDebug("BeamSpotAnalyzer") <<
00104   std::cout <<
00105     "for lumis "<< beginLumiOfBSFit_ << " - " << endLumiOfBSFit_ << std::endl <<
00106     "number of selected tracks = " << theBeamFitter->getNTracks() << std::endl;
00107   std::cout << "number of selected PVs = " << theBeamFitter->getNPVs() << std::endl;
00108   //std::cout << "number of selected PVs per bx: " << std::endl;
00109   //theBeamFitter->getNPVsperBX();
00110 
00111     const edm::TimeValue_t fendtimestamp = lumiSeg.endTime().value();
00112     const std::time_t fendtime = fendtimestamp >> 32;
00113     refBStime[1] = fendtime;
00114 
00115     endLumiOfBSFit_ = lumiSeg.luminosityBlock();
00116     previousLumi_ = endLumiOfBSFit_;
00117 
00118     if ( fitNLumi_ == -1 && resetFitNLumi_ == -1 ) return;
00119 
00120     if (fitNLumi_ > 0 && countLumi_%fitNLumi_!=0) return;
00121 
00122     theBeamFitter->setFitLSRange(beginLumiOfBSFit_,endLumiOfBSFit_);
00123     theBeamFitter->setRefTime(refBStime[0],refBStime[1]);
00124     theBeamFitter->setRun(ftmprun0);
00125 
00126     std::pair<int,int> LSRange = theBeamFitter->getFitLSRange();
00127 
00128     if (theBeamFitter->runPVandTrkFitter()){
00129       reco::BeamSpot bs = theBeamFitter->getBeamSpot();
00130       std::cout << "\n RESULTS OF DEFAULT FIT " << std::endl;
00131       std::cout << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl;
00132       std::cout << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl;
00133       std::cout << " lumi counter # " << countLumi_ << std::endl;
00134       std::cout << bs << std::endl;
00135       std::cout << "[BeamFitter] fit done. \n" << std::endl;
00136     }
00137     else { // Fill in empty beam spot if beamfit fails
00138       reco::BeamSpot bs;
00139       bs.setType(reco::BeamSpot::Fake);
00140       std::cout << "\n Empty Beam spot fit" << std::endl;
00141       std::cout << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl;
00142       std::cout << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl;
00143       std::cout << " lumi counter # " << countLumi_ << std::endl;
00144       std::cout << bs << std::endl;
00145       std::cout << "[BeamFitter] fit failed \n" << std::endl;
00146       //accumulate more events
00147       // dissable this for the moment
00148       //resetFitNLumi_ += 1;
00149       //std::cout << "reset fitNLumi " << resetFitNLumi_ << std::endl;
00150     }
00151 
00152     if (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) {
00153       std::vector<BSTrkParameters> theBSvector = theBeamFitter->getBSvector();
00154       std::cout << "Total number of tracks accumulated = " << theBSvector.size() << std::endl;
00155       std::cout << "Reset track collection for beam fit" <<std::endl;
00156       theBeamFitter->resetTrkVector();
00157       theBeamFitter->resetLSRange();
00158       theBeamFitter->resetCutFlow();
00159       theBeamFitter->resetRefTime();
00160       theBeamFitter->resetPVFitter();
00161       countLumi_=0;
00162       // reset counter to orginal
00163       resetFitNLumi_ = Org_resetFitNLumi_;
00164     }
00165 
00166 }
00167 
00168 
00169 void
00170 BeamSpotAnalyzer::endJob() {
00171   std::cout << "\n-------------------------------------\n" << std::endl;
00172   std::cout << "\n Total number of events processed: "<< ftotalevents << std::endl;
00173   std::cout << "\n-------------------------------------\n\n" << std::endl;
00174 
00175   if ( fitNLumi_ == -1 && resetFitNLumi_ == -1 ) {
00176 
00177           if(theBeamFitter->runPVandTrkFitter()){
00178                   reco::BeamSpot beam_default = theBeamFitter->getBeamSpot();
00179                   std::pair<int,int> LSRange = theBeamFitter->getFitLSRange();
00180 
00181                   std::cout << "\n RESULTS OF DEFAULT FIT:" << std::endl;
00182                   std::cout << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl;
00183                   std::cout << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl;
00184                   std::cout << " lumi counter # " << countLumi_ << std::endl;
00185                   std::cout << beam_default << std::endl;
00186 
00187                   if (write2DB_) {
00188                           std::cout << "\n-------------------------------------\n\n" << std::endl;
00189                           std::cout << " write results to DB..." << std::endl;
00190                           theBeamFitter->write2DB();
00191                   }
00192 
00193                   if (runallfitters_) {
00194                           theBeamFitter->runAllFitter();
00195                   }
00196 
00197           }
00198       if ((runbeamwidthfit_)){
00199                  theBeamFitter->runBeamWidthFitter();
00200                  reco::BeamSpot beam_width = theBeamFitter->getBeamWidth();
00201                  std::cout <<beam_width<< std::endl;
00202       }
00203 
00204           else std::cout << "[BeamSpotAnalyzer] beamfit fails !!!" << std::endl;
00205   }
00206 
00207   std::cout << "[BeamSpotAnalyzer] endJob done \n" << std::endl;
00208 }
00209 
00210 //define this as a plug-in
00211 DEFINE_FWK_MODULE(BeamSpotAnalyzer);