File indexing completed on 2024-04-06 11:59:26
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #include <string>
0015
0016 #include "Calibration/TkAlCaRecoProducers/interface/AlcaBeamSpotProducer.h"
0017 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0018 #include "RecoVertex/BeamSpotProducer/interface/BSFitter.h"
0019
0020 #include "FWCore/Framework/interface/ConsumesCollector.h"
0021 #include "FWCore/Framework/interface/MakerMacros.h"
0022 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0023
0024 #include "CondFormats/BeamSpotObjects/interface/BeamSpotObjects.h"
0025 #include "CondFormats/DataRecord/interface/BeamSpotObjectsRcd.h"
0026 #include "FWCore/Framework/interface/EventSetup.h"
0027 #include "FWCore/Framework/interface/LuminosityBlock.h"
0028
0029 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0030 #include "TMath.h"
0031
0032
0033 AlcaBeamSpotProducer::AlcaBeamSpotProducer(const edm::ParameterSet &iConfig) {
0034
0035 write2DB_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("WriteToDB");
0036 runallfitters_ =
0037 iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunAllFitters");
0038 fitNLumi_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters")
0039 .getUntrackedParameter<int>("fitEveryNLumi", -1);
0040 resetFitNLumi_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters")
0041 .getUntrackedParameter<int>("resetEveryNLumi", -1);
0042 runbeamwidthfit_ =
0043 iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunBeamWidthFit");
0044
0045 theBeamFitter = new BeamFitter(iConfig, consumesCollector());
0046 theBeamFitter->resetTrkVector();
0047 theBeamFitter->resetLSRange();
0048 theBeamFitter->resetCutFlow();
0049 theBeamFitter->resetRefTime();
0050 theBeamFitter->resetPVFitter();
0051
0052 ftotalevents = 0;
0053 ftmprun0 = ftmprun = -1;
0054 countLumi_ = 0;
0055 beginLumiOfBSFit_ = endLumiOfBSFit_ = -1;
0056
0057 produces<reco::BeamSpot, edm::Transition::EndLuminosityBlock>("alcaBeamSpot");
0058 }
0059
0060
0061 AlcaBeamSpotProducer::~AlcaBeamSpotProducer() { delete theBeamFitter; }
0062
0063
0064 void AlcaBeamSpotProducer::produce(edm::Event &iEvent, const edm::EventSetup &iSetup) {
0065 ftotalevents++;
0066 theBeamFitter->readEvent(iEvent);
0067 ftmprun = iEvent.id().run();
0068 }
0069
0070
0071 void AlcaBeamSpotProducer::beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) {
0072 const edm::TimeValue_t fbegintimestamp = lumiSeg.beginTime().value();
0073 const std::time_t ftmptime = fbegintimestamp >> 32;
0074
0075 if (countLumi_ == 0 || (resetFitNLumi_ > 0 && countLumi_ % resetFitNLumi_ == 0)) {
0076 ftmprun0 = lumiSeg.run();
0077 ftmprun = ftmprun0;
0078 beginLumiOfBSFit_ = lumiSeg.luminosityBlock();
0079 refBStime[0] = ftmptime;
0080 }
0081
0082 countLumi_++;
0083 }
0084
0085
0086 void AlcaBeamSpotProducer::endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) {}
0087
0088
0089 void AlcaBeamSpotProducer::endLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) {
0090 const edm::TimeValue_t fendtimestamp = lumiSeg.endTime().value();
0091 const std::time_t fendtime = fendtimestamp >> 32;
0092 refBStime[1] = fendtime;
0093
0094 endLumiOfBSFit_ = lumiSeg.luminosityBlock();
0095
0096 if (fitNLumi_ == -1 && resetFitNLumi_ == -1)
0097 return;
0098
0099 if (fitNLumi_ > 0 && countLumi_ % fitNLumi_ != 0)
0100 return;
0101
0102 theBeamFitter->setFitLSRange(beginLumiOfBSFit_, endLumiOfBSFit_);
0103 theBeamFitter->setRefTime(refBStime[0], refBStime[1]);
0104 theBeamFitter->setRun(ftmprun0);
0105
0106 std::pair<int, int> LSRange = theBeamFitter->getFitLSRange();
0107
0108 reco::BeamSpot bs;
0109 if (theBeamFitter->runPVandTrkFitter()) {
0110 bs = theBeamFitter->getBeamSpot();
0111 edm::LogInfo("AlcaBeamSpotProducer") << "\n RESULTS OF DEFAULT FIT " << std::endl
0112 << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
0113 << " for lumi blocks : " << LSRange.first << " - " << LSRange.second
0114 << std::endl
0115 << " lumi counter # " << countLumi_ << std::endl
0116 << bs << std::endl
0117 << "fit done. \n"
0118 << std::endl;
0119 } else {
0120 bs.setType(reco::BeamSpot::Fake);
0121 edm::LogInfo("AlcaBeamSpotProducer") << "\n Empty Beam spot fit" << std::endl
0122 << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
0123 << " for lumi blocks : " << LSRange.first << " - " << LSRange.second
0124 << std::endl
0125 << " lumi counter # " << countLumi_ << std::endl
0126 << bs << std::endl
0127 << "fit failed \n"
0128 << std::endl;
0129 }
0130
0131 auto result = std::make_unique<reco::BeamSpot>();
0132 *result = bs;
0133 lumiSeg.put(std::move(result), std::string("alcaBeamSpot"));
0134
0135 if (resetFitNLumi_ > 0 && countLumi_ % resetFitNLumi_ == 0) {
0136 std::vector<BSTrkParameters> theBSvector = theBeamFitter->getBSvector();
0137 edm::LogInfo("AlcaBeamSpotProducer") << "Total number of tracks accumulated = " << theBSvector.size() << std::endl
0138 << "Reset track collection for beam fit" << std::endl;
0139 theBeamFitter->resetTrkVector();
0140 theBeamFitter->resetLSRange();
0141 theBeamFitter->resetCutFlow();
0142 theBeamFitter->resetRefTime();
0143 theBeamFitter->resetPVFitter();
0144 countLumi_ = 0;
0145 }
0146 }
0147
0148 DEFINE_FWK_MODULE(AlcaBeamSpotProducer);