CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlcaBeamSpotProducer.cc
Go to the documentation of this file.
1 
15 // C++ standard
16 #include <string>
17 // CMS
21 
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);
49 
50  ftotalevents = 0;
51  ftmprun0 = ftmprun = -1;
52  countLumi_ = 0;
54 
55  produces<reco::BeamSpot, edm::InLumi>("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  const edm::TimeValue_t fendtimestamp = lumiSeg.endTime().value();
89  const std::time_t fendtime = fendtimestamp >> 32;
90  refBStime[1] = fendtime;
91 
92  endLumiOfBSFit_ = lumiSeg.luminosityBlock();
93 
94  if ( fitNLumi_ == -1 && resetFitNLumi_ == -1 ) return;
95 
96  if (fitNLumi_ > 0 && countLumi_%fitNLumi_!=0) return;
97 
101 
102  int * LSRange = theBeamFitter->getFitLSRange();
103 
104  reco::BeamSpot bs;
106  bs = theBeamFitter->getBeamSpot();
107  edm::LogInfo("AlcaBeamSpotProducer")
108  << "\n RESULTS OF DEFAULT FIT " << std::endl
109  << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
110  << " for lumi blocks : " << LSRange[0] << " - " << LSRange[1] << std::endl
111  << " lumi counter # " << countLumi_ << std::endl
112  << bs << std::endl
113  << "fit done. \n" << std::endl;
114  }
115  else { // Fill in empty beam spot if beamfit fails
117  edm::LogInfo("AlcaBeamSpotProducer")
118  << "\n Empty Beam spot fit" << std::endl
119  << " for runs: " << ftmprun0 << " - " << ftmprun << std::endl
120  << " for lumi blocks : " << LSRange[0] << " - " << LSRange[1] << std::endl
121  << " lumi counter # " << countLumi_ << std::endl
122  << bs << std::endl
123  << "fit failed \n" << std::endl;
124  }
125 
126  std::auto_ptr<reco::BeamSpot> result(new reco::BeamSpot);
127  *result = bs;
128  lumiSeg.put(result, std::string("alcaBeamSpot"));
129 
130  if (resetFitNLumi_ > 0 && countLumi_%resetFitNLumi_ == 0) {
131  std::vector<BSTrkParameters> theBSvector = theBeamFitter->getBSvector();
132  edm::LogInfo("AlcaBeamSpotProducer")
133  << "Total number of tracks accumulated = " << theBSvector.size() << std::endl
134  << "Reset track collection for beam fit" <<std::endl;
140  countLumi_=0;
141  }
142 }
143 
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
bool runPVandTrkFitter()
Definition: BeamFitter.cc:394
void setRun(int run)
Definition: BeamFitter.h:80
std::vector< BSTrkParameters > getBSvector()
Definition: BeamFitter.h:62
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void setRefTime(std::time_t t0, std::time_t t1)
Definition: BeamFitter.h:52
void resetTrkVector()
Definition: BeamFitter.h:48
Timestamp const & beginTime() const
void readEvent(const edm::Event &iEvent)
Definition: BeamFitter.cc:211
void setType(BeamType type)
set beam type
Definition: BeamSpot.h:122
LuminosityBlockNumber_t luminosityBlock() const
int iEvent
Definition: GenABIO.cc:243
AlcaBeamSpotProducer(const edm::ParameterSet &)
virtual void produce(edm::Event &iEvent, const edm::EventSetup &iSetup)
Timestamp const & endTime() const
tuple result
Definition: query.py:137
RunNumber_t run() const
int * getFitLSRange()
Definition: BeamFitter.h:70
void setFitLSRange(int ls0, int ls1)
Definition: BeamFitter.h:76
void resetLSRange()
Definition: BeamFitter.h:50
unsigned long long TimeValue_t
Definition: Timestamp.h:27
TimeValue_t value() const
Definition: Timestamp.cc:72
edm::EventID id() const
Definition: EventBase.h:56
void resetPVFitter()
Definition: BeamFitter.h:56
void resetCutFlow()
Definition: BeamFitter.h:64
virtual void beginLuminosityBlock(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup)
virtual void endLuminosityBlock(edm::LuminosityBlock &lumiSeg, const edm::EventSetup &iSetup)
reco::BeamSpot getBeamSpot()
Definition: BeamFitter.h:60
void put(std::auto_ptr< PROD > product)
Put a new product.
void resetRefTime()
Definition: BeamFitter.h:51