CMS 3D CMS Logo

AlcaBeamSpotProducer.cc
Go to the documentation of this file.
1 
14 // C++ standard
15 #include <string>
16 // CMS
20 
24 
30 
32 #include "TMath.h"
33 
34 //--------------------------------------------------------------------------------------------------
36  // get parameter
37  write2DB_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("WriteToDB");
38  runallfitters_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunAllFitters");
39  fitNLumi_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getUntrackedParameter<int>("fitEveryNLumi",-1);
40  resetFitNLumi_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getUntrackedParameter<int>("resetEveryNLumi",-1);
41  runbeamwidthfit_ = iConfig.getParameter<edm::ParameterSet>("AlcaBeamSpotProducerParameters").getParameter<bool>("RunBeamWidthFit");
42 
43  theBeamFitter = new BeamFitter(iConfig, consumesCollector());
49 
50  ftotalevents = 0;
51  ftmprun0 = ftmprun = -1;
52  countLumi_ = 0;
54 
55  produces<reco::BeamSpot, edm::Transition::EndLuminosityBlock>("alcaBeamSpot");
56 }
57 
58 //--------------------------------------------------------------------------------------------------
60  delete theBeamFitter;
61 }
62 
63 //--------------------------------------------------------------------------------------------------
65  ftotalevents++;
66  theBeamFitter->readEvent(iEvent);
67  ftmprun = iEvent.id().run();
68 }
69 
70 //--------------------------------------------------------------------------------------------------
72  const edm::TimeValue_t fbegintimestamp = lumiSeg.beginTime().value();
73  const std::time_t ftmptime = fbegintimestamp >> 32;
74 
75  if ( countLumi_ == 0 || (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) ) {
76  ftmprun0 = lumiSeg.run();
77  ftmprun = ftmprun0;
79  refBStime[0] = ftmptime;
80  }
81 
82  countLumi_++;
83 
84 }
85 
86 //--------------------------------------------------------------------------------------------------
88 }
89 
90 //--------------------------------------------------------------------------------------------------
92  const edm::TimeValue_t fendtimestamp = lumiSeg.endTime().value();
93  const std::time_t fendtime = fendtimestamp >> 32;
94  refBStime[1] = fendtime;
95 
96  endLumiOfBSFit_ = lumiSeg.luminosityBlock();
97 
98  if ( fitNLumi_ == -1 && resetFitNLumi_ == -1 ) return;
99 
100  if (fitNLumi_ > 0 && countLumi_%fitNLumi_!=0) return;
101 
105 
106  std::pair<int,int> LSRange = theBeamFitter->getFitLSRange();
107 
108  reco::BeamSpot bs;
110  bs = theBeamFitter->getBeamSpot();
111  edm::LogInfo("AlcaBeamSpotProducer")
112  << "\n RESULTS OF DEFAULT FIT " << std::endl
113  << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
114  << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl
115  << " lumi counter # " << countLumi_ << std::endl
116  << bs << std::endl
117  << "fit done. \n" << std::endl;
118  }
119  else { // Fill in empty beam spot if beamfit fails
121  edm::LogInfo("AlcaBeamSpotProducer")
122  << "\n Empty Beam spot fit" << std::endl
123  << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
124  << " for lumi blocks : " << LSRange.first << " - " << LSRange.second << std::endl
125  << " lumi counter # " << countLumi_ << std::endl
126  << bs << std::endl
127  << "fit failed \n" << std::endl;
128  }
129 
130  auto result = std::make_unique<reco::BeamSpot>();
131  *result = bs;
132  lumiSeg.put(std::move(result), std::string("alcaBeamSpot"));
133 
134  if (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) {
135  std::vector<BSTrkParameters> theBSvector = theBeamFitter->getBSvector();
136  edm::LogInfo("AlcaBeamSpotProducer")
137  << "Total number of tracks accumulated = " << theBSvector.size() << std::endl
138  << "Reset track collection for beam fit" <<std::endl;
144  countLumi_=0;
145  }
146 }
147 
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
bool runPVandTrkFitter()
Definition: BeamFitter.cc:398
void setRun(int run)
Definition: BeamFitter.h:124
std::vector< BSTrkParameters > getBSvector()
Definition: BeamFitter.h:96
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void resetTrkVector()
Definition: BeamFitter.h:53
Timestamp const & beginTime() const
void readEvent(const edm::Event &iEvent)
Definition: BeamFitter.cc:216
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:131
LuminosityBlockNumber_t luminosityBlock() const
int iEvent
Definition: GenABIO.cc:230
AlcaBeamSpotProducer(const edm::ParameterSet &)
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
void put(std::unique_ptr< PROD > product)
Put a new product.
Timestamp const & endTime() const
RunNumber_t run() const
void setFitLSRange(int ls0, int ls1)
Definition: BeamFitter.h:120
void resetLSRange()
Definition: BeamFitter.h:55
unsigned long long TimeValue_t
Definition: Timestamp.h:28
void endLuminosityBlockProduce(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup) final
edm::EventID id() const
Definition: EventBase.h:60
void resetPVFitter()
Definition: BeamFitter.h:69
void produce(edm::Event &iEvent, const edm::EventSetup &iSetup) final
void resetCutFlow()
Definition: BeamFitter.h:105
void setRefTime(time_t t0, time_t t1)
Definition: BeamFitter.h:57
std::pair< int, int > getFitLSRange()
Definition: BeamFitter.h:117
TimeValue_t value() const
Definition: Timestamp.h:56
reco::BeamSpot getBeamSpot()
Definition: BeamFitter.h:94
def move(src, dest)
Definition: eostools.py:510
void endLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
void beginLuminosityBlock(edm::LuminosityBlock const &lumiSeg, const edm::EventSetup &iSetup) final
void resetRefTime()
Definition: BeamFitter.h:56