CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_8_patch3/src/Calibration/TkAlCaRecoProducers/plugins/AlcaBeamSpotProducer.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 "Calibration/TkAlCaRecoProducers/interface/AlcaBeamSpotProducer.h"
00020 #include "RecoVertex/BeamSpotProducer/interface/BSFitter.h"
00021 
00022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00023 #include "FWCore/Framework/interface/MakerMacros.h"
00024 
00025 #include "FWCore/Framework/interface/ESHandle.h"
00026 #include "FWCore/Framework/interface/EventSetup.h"
00027 #include "FWCore/Framework/interface/LuminosityBlock.h"
00028 #include "CondFormats/DataRecord/interface/BeamSpotObjectsRcd.h"
00029 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
00030 
00031 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00032 #include "TMath.h"
00033 
00034 //--------------------------------------------------------------------------------------------------
00035 AlcaBeamSpotProducer::AlcaBeamSpotProducer(const edm::ParameterSet& iConfig){
00036   // get parameter
00037   write2DB_        = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("WriteToDB");
00038   runallfitters_   = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunAllFitters");
00039   fitNLumi_        = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getUntrackedParameter<int>("fitEveryNLumi",-1);
00040   resetFitNLumi_   = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getUntrackedParameter<int>("resetEveryNLumi",-1);
00041   runbeamwidthfit_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunBeamWidthFit");
00042   
00043   theBeamFitter = new BeamFitter(iConfig);
00044   theBeamFitter->resetTrkVector();
00045   theBeamFitter->resetLSRange();
00046   theBeamFitter->resetCutFlow();
00047   theBeamFitter->resetRefTime();
00048   theBeamFitter->resetPVFitter();
00049 
00050   ftotalevents = 0;
00051   ftmprun0 = ftmprun = -1;
00052   countLumi_ = 0;
00053   beginLumiOfBSFit_ = endLumiOfBSFit_ = -1;
00054   
00055   produces<reco::BeamSpot, edm::InLumi>("alcaBeamSpot");
00056 }
00057 
00058 //--------------------------------------------------------------------------------------------------
00059 AlcaBeamSpotProducer::~AlcaBeamSpotProducer(){
00060   delete theBeamFitter;
00061 }
00062 
00063 //--------------------------------------------------------------------------------------------------
00064 void AlcaBeamSpotProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup){
00065   ftotalevents++;
00066   theBeamFitter->readEvent(iEvent);
00067   ftmprun = iEvent.id().run();
00068 }
00069 
00070 //--------------------------------------------------------------------------------------------------
00071 void AlcaBeamSpotProducer::beginLuminosityBlock(edm::LuminosityBlock& lumiSeg, const edm::EventSetup& iSetup){
00072   const edm::TimeValue_t fbegintimestamp = lumiSeg.beginTime().value();
00073   const std::time_t ftmptime = fbegintimestamp >> 32;
00074 
00075   if ( countLumi_ == 0 || (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) ) {
00076     ftmprun0 = lumiSeg.run();
00077     ftmprun = ftmprun0;
00078     beginLumiOfBSFit_ = lumiSeg.luminosityBlock();
00079     refBStime[0] = ftmptime;
00080   }
00081     
00082   countLumi_++;
00083   
00084 }
00085 
00086 //--------------------------------------------------------------------------------------------------
00087 void AlcaBeamSpotProducer::endLuminosityBlock(edm::LuminosityBlock& lumiSeg, const edm::EventSetup& iSetup){
00088   const edm::TimeValue_t fendtimestamp = lumiSeg.endTime().value();
00089   const std::time_t fendtime = fendtimestamp >> 32;
00090   refBStime[1] = fendtime;
00091     
00092   endLumiOfBSFit_ = lumiSeg.luminosityBlock();
00093     
00094   if ( fitNLumi_ == -1 && resetFitNLumi_ == -1 ) return;
00095         
00096   if (fitNLumi_ > 0 && countLumi_%fitNLumi_!=0) return;
00097 
00098   theBeamFitter->setFitLSRange(beginLumiOfBSFit_,endLumiOfBSFit_);
00099   theBeamFitter->setRefTime(refBStime[0],refBStime[1]);
00100   theBeamFitter->setRun(ftmprun0);
00101     
00102   std::pair<int,int> LSRange = theBeamFitter->getFitLSRange();
00103 
00104   reco::BeamSpot bs;
00105   if (theBeamFitter->runPVandTrkFitter()){
00106     bs = theBeamFitter->getBeamSpot();
00107     edm::LogInfo("AlcaBeamSpotProducer")
00108         << "\n RESULTS OF DEFAULT FIT " << std::endl
00109         << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
00110         << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl
00111         << " lumi counter # " << countLumi_ << std::endl
00112         << bs << std::endl
00113         << "fit done. \n" << std::endl; 
00114   }
00115   else { // Fill in empty beam spot if beamfit fails
00116     bs.setType(reco::BeamSpot::Fake);
00117     edm::LogInfo("AlcaBeamSpotProducer")
00118         << "\n Empty Beam spot fit" << std::endl
00119         << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
00120         << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl
00121         << " lumi counter # " << countLumi_ << std::endl
00122         << bs << std::endl
00123         << "fit failed \n" << std::endl;
00124   }
00125 
00126   std::auto_ptr<reco::BeamSpot> result(new reco::BeamSpot);
00127   *result = bs;
00128   lumiSeg.put(result, std::string("alcaBeamSpot"));
00129         
00130   if (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) {
00131     std::vector<BSTrkParameters> theBSvector = theBeamFitter->getBSvector();
00132     edm::LogInfo("AlcaBeamSpotProducer")
00133         << "Total number of tracks accumulated = " << theBSvector.size() << std::endl
00134         << "Reset track collection for beam fit" <<std::endl;
00135     theBeamFitter->resetTrkVector();
00136     theBeamFitter->resetLSRange();
00137     theBeamFitter->resetCutFlow();
00138     theBeamFitter->resetRefTime();
00139     theBeamFitter->resetPVFitter();
00140     countLumi_=0;
00141   }
00142 }
00143 
00144 DEFINE_FWK_MODULE(AlcaBeamSpotProducer);