CMS 3D CMS Logo

AlcaBeamSpotManager.cc
Go to the documentation of this file.
1 
11 #include <vector>
12 #include <cmath>
13 #include <climits>
14 
15 using namespace edm;
16 using namespace reco;
17 using namespace std;
18 
19 //--------------------------------------------------------------------------------------------------
21 }
22 
23 //--------------------------------------------------------------------------------------------------
25  beamSpotOutputBase_(iConfig.getParameter<ParameterSet>("AlcaBeamSpotHarvesterParameters").getUntrackedParameter<std::string>("BeamSpotOutputBase")),
26  beamSpotModuleName_(iConfig.getParameter<ParameterSet>("AlcaBeamSpotHarvesterParameters").getUntrackedParameter<std::string>("BeamSpotModuleName")),
27  beamSpotLabel_ (iConfig.getParameter<ParameterSet>("AlcaBeamSpotHarvesterParameters").getUntrackedParameter<std::string>("BeamSpotLabel")),
28  sigmaZCut_ (iConfig.getParameter<ParameterSet>("AlcaBeamSpotHarvesterParameters").getUntrackedParameter<double>("SigmaZCut"))
29 {
32  LogInfo("AlcaBeamSpotManager")
33  << "Output base: " << beamSpotOutputBase_
34  << std::endl;
35  reset();
36 }
37 
38 //--------------------------------------------------------------------------------------------------
40 }
41 
42 //--------------------------------------------------------------------------------------------------
44  beamSpotMap_.clear();
45 }
46 //--------------------------------------------------------------------------------------------------
48 
49  Handle<BeamSpot> beamSpotHandle;
50  iLumi.getByToken(beamSpotToken_, beamSpotHandle);
51 
52  if(beamSpotHandle.isValid()) { // check the product
53  std::pair<edm::Timestamp,reco::BeamSpot> time_bs (iLumi.beginTime(), *beamSpotHandle);
54  beamSpotMap_[iLumi.luminosityBlock()] = time_bs;
55  const BeamSpot* aBeamSpot = &beamSpotMap_[iLumi.luminosityBlock()].second;
56  aBeamSpot = beamSpotHandle.product();
57  LogInfo("AlcaBeamSpotManager")
58  << "Lumi: " << iLumi.luminosityBlock() << std::endl;
59  LogInfo("AlcaBeamSpotManager")
60  << *aBeamSpot << std::endl;
61  }
62  else {
63  LogInfo("AlcaBeamSpotManager")
64  << "Lumi: " << iLumi.luminosityBlock() << std::endl;
65  LogInfo("AlcaBeamSpotManager")
66  << " BS is not valid!" << std::endl;
67  }
68 
69 }
70 
71 //--------------------------------------------------------------------------------------------------
73  vector<bsMap_iterator> listToErase;
74  for(bsMap_iterator it=beamSpotMap_.begin(); it!=beamSpotMap_.end();it++){
75  if(it->second.second.type() != BeamSpot::Tracker || it->second.second.sigmaZ()<sigmaZCut_ ) {
76  listToErase.push_back(it);
77  }
78  }
79  for(vector<bsMap_iterator>::iterator it=listToErase.begin(); it !=listToErase.end(); it++){
80  beamSpotMap_.erase(*it);
81  }
82  if(beamSpotMap_.size() <= 1){
83  return;
84  }
85  //Return only if lumibased since the collapsing alghorithm requires the next and next to next lumi sections
86  else if(beamSpotMap_.size() == 2 && beamSpotOutputBase_ == "lumibased"){
87  return;
88  }
89  if(beamSpotOutputBase_ == "lumibased"){
90 // bsMap_iterator referenceBS = beamSpotMap_.begin();
91  bsMap_iterator firstBS = beamSpotMap_.begin();
92 // bsMap_iterator lastBS = beamSpotMap_.begin();
93  bsMap_iterator currentBS = beamSpotMap_.begin();
94  bsMap_iterator nextBS = ++beamSpotMap_.begin();
95  bsMap_iterator nextNextBS = ++(++(beamSpotMap_.begin()));
96 
97  reco::BeamSpot currentBSObj;
98  reco::BeamSpot nextBSObj;
99 
100  map<LuminosityBlockNumber_t,std::pair<edm::Timestamp, BeamSpot>> tmpBeamSpotMap_;
101  bool docreate = true;
102  bool endOfRun = false;//Added
103  bool docheck = true;//Added
104  bool foundShift = false;
105  long countlumi = 0;//Added
106  string tmprun = "";//Added
107  long maxNlumis = 20;//Added
108 // if weighted:
109 // maxNlumis = 999999999
110 
111 
112  unsigned int iteration = 0;
113  //while(nextNextBS!=beamSpotMap_.end()){
114  while(nextBS!=beamSpotMap_.end()){
115  LogInfo("AlcaBeamSpotManager")
116  << "Iteration: " << iteration << " size: " << beamSpotMap_.size() << "\n"
117  << "Lumi: " << currentBS->first << "\n"
118  << currentBS->second.second
119  << "\n" << nextBS->first << "\n" << nextBS->second.second
120  << endl;
121  if (nextNextBS!=beamSpotMap_.end())
122  LogInfo("AlcaBeamSpotManager")
123  << nextNextBS->first << "\n" << nextNextBS->second.second
124  << endl;
125 
126 
127  if(docreate){
128  firstBS = currentBS;
129  docreate = false;//Added
130  }
131  //if(iteration >= beamSpotMap_.size()-3){
132  if(iteration >= beamSpotMap_.size()-2){
133  LogInfo("AlcaBeamSpotManager")
134  << "Reached lumi " << currentBS->first
135  << " now close payload because end of data has been reached.";
136  docreate = true;
137  endOfRun = true;
138  }
139  // check we run over the same run
140 // if (ibeam->first.Run() != inextbeam->first.Run()){
141 // LogInfo("AlcaBeamSpotManager")
142 // << "close payload because end of run.";
143 // docreate = true;
144 // }
145  // check maximum lumi counts
146  if (countlumi == maxNlumis -1){
147  LogInfo("AlcaBeamSpotManager")
148  << "close payload because maximum lumi sections accumulated within run ";
149  docreate = true;
150  countlumi = 0;
151  }
152 // if weighted:
153 // docheck = False
154  // check offsets
155  if(docheck){
156  foundShift = false;
157  LogInfo("AlcaBeamSpotManager")
158  << "Checking checking!" << endl;
159  float limit = 0;
160  pair<float,float> adelta1;
161  pair<float,float> adelta2;
162  pair<float,float> adelta;
163  pair<float,float> adelta1dxdz;
164  pair<float,float> adelta2dxdz;
165  pair<float,float> adelta1dydz;
166  pair<float,float> adelta2dydz;
167  pair<float,float> adelta1widthX;
168  pair<float,float> adelta2widthX;
169  pair<float,float> adelta1widthY;
170  pair<float,float> adelta2widthY;
171  pair<float,float> adelta1z0;
172  pair<float,float> adelta1sigmaZ;
173 
174  // define minimum limit
175  float min_limit = 0.0025;
176 
177  // limit for x and y
178  limit = currentBS->second.second.BeamWidthX()/2.;
179  if(limit < min_limit){
180  limit = min_limit;
181  }
182 
183  currentBSObj = currentBS->second.second;
184  nextBSObj = nextBS->second.second;
185 
186  //check movements in X
187  adelta1 = delta(currentBSObj.x0(), currentBSObj.x0Error(), nextBSObj.x0(), nextBSObj.x0Error());
188  adelta2 = pair<float,float>(0.,1.e9);
189  if (nextNextBS->second.second.type() != -1){
190  adelta2 = delta(nextBS->second.second.x0(), nextBSObj.x0Error(), nextNextBS->second.second.x0(), nextNextBS->second.second.x0Error());
191  }
192  bool deltaX = (deltaSig(adelta1.first,adelta1.second) > 3.5 && adelta1.first >= limit)?true:false;
193  if(iteration < beamSpotMap_.size()-2){
194  if( !deltaX && adelta1.first*adelta2.first > 0. && fabs(adelta1.first+adelta2.first) >= limit){
195  LogInfo("AlcaBeamSpotManager")
196  << " positive, " << (adelta1.first+adelta2.first) << " limit=" << limit << endl;
197  deltaX = true;
198  }
199  else if( deltaX && adelta1.first*adelta2.first < 0 && adelta2.first != 0 && fabs(adelta1.first/adelta2.first) > 0.33 && fabs(adelta1.first/adelta2.first) < 3){
200  LogInfo("AlcaBeamSpotManager")
201  << " negative, " << adelta1.first/adelta2.first << endl;
202  deltaX = false;
203  }
204  }
205 
206  //calculating all deltas
207  adelta1dxdz = delta(currentBSObj.dxdz(), currentBSObj.dxdzError(), nextBSObj.dxdz(), nextBSObj.dxdzError());
208  adelta2dxdz = pair<float,float>(0.,1.e9);
209  adelta1dydz = delta(currentBSObj.dydz(), currentBSObj.dydzError(), nextBSObj.dydz(), nextBSObj.dydzError());
210  adelta2dydz = pair<float,float>(0.,1.e9);
211  adelta1widthX = delta(currentBSObj.BeamWidthX(), currentBSObj.BeamWidthXError(), nextBSObj.BeamWidthX(), nextBSObj.BeamWidthXError());
212  adelta2widthX = pair<float,float>(0.,1.e9);
213  adelta1widthY = delta(currentBSObj.BeamWidthY(), currentBSObj.BeamWidthYError(), nextBSObj.BeamWidthY(), nextBSObj.BeamWidthYError());
214  adelta2widthY = pair<float,float>(0.,1.e9);
215  adelta1z0 = delta(currentBSObj.z0(), currentBSObj.z0Error(), nextBSObj.z0(), nextBSObj.z0Error());
216  adelta1sigmaZ = delta(currentBSObj.sigmaZ(), currentBSObj.sigmaZ0Error(), nextBSObj.sigmaZ(), nextBSObj.sigmaZ0Error());
217 
218  //check movements in Y
219  adelta1 = delta(currentBSObj.y0(), currentBSObj.y0Error(), nextBSObj.y0(), nextBSObj.y0Error());
220  adelta2 = pair<float,float>(0.,1.e9);
221  if( nextNextBS->second.second.type() != BeamSpot::Unknown){
222  adelta2 = delta(nextBSObj.y0(), nextBSObj.y0Error(), nextNextBS->second.second.y0(), nextNextBS->second.second.y0Error());
223  adelta2dxdz = delta(nextBSObj.dxdz(), nextBSObj.dxdzError(), nextNextBS->second.second.dxdz(), nextNextBS->second.second.dxdzError());
224  adelta2dydz = delta(nextBSObj.dydz(), nextBSObj.dydzError(), nextNextBS->second.second.dydz(), nextNextBS->second.second.dydzError());
225  adelta2widthX = delta(nextBSObj.BeamWidthX(), nextBSObj.BeamWidthXError(), nextNextBS->second.second.BeamWidthX(), nextNextBS->second.second.BeamWidthXError());
226  adelta2widthY = delta(nextBSObj.BeamWidthY(), nextBSObj.BeamWidthYError(), nextNextBS->second.second.BeamWidthY(), nextNextBS->second.second.BeamWidthYError());
227 
228  }
229  bool deltaY = (deltaSig(adelta1.first,adelta1.second) > 3.5 && adelta1.first >= limit)?true:false;
230  if(iteration < beamSpotMap_.size()-2){
231  if( !deltaY && adelta1.first*adelta2.first > 0. && fabs(adelta1.first+adelta2.first) >= limit){
232  LogInfo("AlcaBeamSpotManager")
233  << " positive, " << (adelta1.first+adelta2.first) << " limit=" << limit << endl;
234  deltaY = true;
235  }
236  else if( deltaY && adelta1.first*adelta2.first < 0 && adelta2.first != 0 && fabs(adelta1.first/adelta2.first) > 0.33 && fabs(adelta1.first/adelta2.first) < 3){
237  LogInfo("AlcaBeamSpotManager")
238  << " negative, " << adelta1.first/adelta2.first << endl;
239  deltaY = false;
240  }
241  }
242 
243  limit = currentBSObj.sigmaZ()/2.;
244  bool deltaZ = (deltaSig(adelta1z0.first,adelta1z0.second) > 3.5 && fabs(adelta1z0.first) >= limit)?true:false;
245  adelta = delta(currentBSObj.sigmaZ(), currentBSObj.sigmaZ0Error(), nextBSObj.sigmaZ(), nextBSObj.sigmaZ0Error());
246  bool deltasigmaZ = (deltaSig(adelta.first,adelta.second) > 5.0)?true:false;
247  bool deltadxdz = false;
248  bool deltadydz = false;
249  bool deltawidthX = false;
250  bool deltawidthY = false;
251 
252  if(iteration < beamSpotMap_.size()-2){
253 
254  adelta = delta(currentBSObj.dxdz(), currentBSObj.dxdzError(), nextBSObj.dxdz(), nextBSObj.dxdzError());
255  deltadxdz = (deltaSig(adelta.first,adelta.second) > 5.0)?true:false;
256  if(deltadxdz && (adelta1dxdz.first*adelta2dxdz.first) < 0 && adelta2dxdz.first != 0 && fabs(adelta1dxdz.first/adelta2dxdz.first) > 0.33 && fabs(adelta1dxdz.first/adelta2dxdz.first) < 3){
257  deltadxdz = false;
258  }
259 
260  adelta = delta(currentBSObj.dydz(), currentBSObj.dydzError(), nextBSObj.dydz(), nextBSObj.dydzError());
261  deltadydz = (deltaSig(adelta.first,adelta.second) > 5.0)?true:false;
262  if(deltadydz && (adelta1dydz.first*adelta2dydz.first) < 0 && adelta2dydz.first != 0 && fabs(adelta1dydz.first/adelta2dydz.first) > 0.33 && fabs(adelta1dydz.first/adelta2dydz.first) < 3){
263  deltadydz = false;
264  }
265 
266  adelta = delta(currentBSObj.BeamWidthX(), currentBSObj.BeamWidthXError(), nextBSObj.BeamWidthX(), nextBSObj.BeamWidthXError());
267  deltawidthX = (deltaSig(adelta.first,adelta.second) > 5.0)?true:false;
268  if(deltawidthX && (adelta1widthX.first*adelta2widthX.first) < 0 && adelta2widthX.first != 0 && fabs(adelta1widthX.first/adelta2widthX.first) > 0.33 && fabs(adelta1widthX.first/adelta2widthX.first) < 3){
269  deltawidthX = false;
270  }
271 
272  adelta = delta(currentBSObj.BeamWidthY(), currentBSObj.BeamWidthYError(), nextBSObj.BeamWidthY(), nextBSObj.BeamWidthYError());
273  deltawidthY = (deltaSig(adelta.first,adelta.second) > 5.0)?true:false;
274  if(deltawidthY && (adelta1widthY.first*adelta2widthY.first) < 0 && adelta2widthY.first != 0 && fabs(adelta1widthY.first/adelta2widthY.first) > 0.33 && fabs(adelta1widthY.first/adelta2widthY.first) < 3){
275  deltawidthY = false;
276  }
277 
278  }
279  if (deltaX || deltaY || deltaZ || deltasigmaZ || deltadxdz || deltadydz || deltawidthX || deltawidthY){
280  docreate = true;
281  foundShift = true;
282  LogInfo("AlcaBeamSpotManager")
283  << "close payload because of movement in"
284  << " X=" << deltaX
285  << ", Y=" << deltaY
286  << ", Z=" << deltaZ
287  << ", sigmaZ=" << deltasigmaZ
288  << ", dxdz=" << deltadxdz
289  << ", dydz=" << deltadydz
290  << ", widthX=" << deltawidthX
291  << ", widthY=" << deltawidthY
292  << endl;
293  }
294 
295  }
296  if(docreate){
297  std::pair<edm::Timestamp,reco::BeamSpot> thepair;
298  if(foundShift){
299  thepair = std::make_pair( currentBS->second.first, weight(firstBS,nextBS));
300  tmpBeamSpotMap_[firstBS->first] = thepair;
301  if (endOfRun){
302  //if we're here, then we need to found a shift in the last LS
303  //We already created a new IOV, now create one just for the last LS
304  thepair = std::make_pair( nextBS->second.first, nextBSObj);
305  tmpBeamSpotMap_[nextBS->first] = thepair;
306  }
307  }
308  else if(!foundShift && !endOfRun){ //maxLS reached
309  thepair = std::make_pair( currentBS->second.first, weight(firstBS,nextBS));
310  tmpBeamSpotMap_[firstBS->first] = thepair;
311  }
312  else { // end of run with no shift detectred in last LS
313  thepair = std::make_pair( currentBS->second.first, weight(firstBS,beamSpotMap_.end()));
314  tmpBeamSpotMap_[firstBS->first] = thepair;
315  }
316  firstBS = nextBS;
317  countlumi = 0;
318 
319  }
320  //tmprun = currentBS->second.Run
321  // increase the counter by one only if the IOV hasn't been closed
322  if (!docreate) ++countlumi;
323 
324  currentBS = nextBS;
325  nextBS = nextNextBS;
326  nextNextBS++;
327  ++iteration;
328  }
329  beamSpotMap_.clear();
330  beamSpotMap_ = tmpBeamSpotMap_;
331  }
332  else if(beamSpotOutputBase_ == "runbased"){
333  LuminosityBlockNumber_t firstLumi = beamSpotMap_.begin()->first;
334  std::pair<edm::Timestamp,reco::BeamSpot> thepair (beamSpotMap_.begin()->second.first, weight(beamSpotMap_.begin(),beamSpotMap_.end()));
335  beamSpotMap_.clear();
336  beamSpotMap_[firstLumi] = thepair;
337  }
338  else{
339  LogInfo("AlcaBeamSpotManager")
340  << "Unrecognized BeamSpotOutputBase parameter: " << beamSpotOutputBase_
341  << endl;
342  }
343 }
344 
345 //--------------------------------------------------------------------------------------------------
347  const bsMap_iterator& end){
348  double x,xError = 0;
349  double y,yError = 0;
350  double z,zError = 0;
351  double sigmaZ,sigmaZError = 0;
352  double dxdz,dxdzError = 0;
353  double dydz,dydzError = 0;
354  double widthX,widthXError = 0;
355  double widthY,widthYError = 0;
356  LogInfo("AlcaBeamSpotManager")
357  << "Weighted BeamSpot will span lumi "
358  << begin->first << " to " << end->first
359  << endl;
360 
362  for(bsMap_iterator it=begin; it!=end; it++){
363  weight(x , xError , it->second.second.x0() , it->second.second.x0Error());
364  weight(y , yError , it->second.second.y0() , it->second.second.y0Error());
365  weight(z , zError , it->second.second.z0() , it->second.second.z0Error());
366  weight(sigmaZ, sigmaZError, it->second.second.sigmaZ() , it->second.second.sigmaZ0Error());
367  weight(dxdz , dxdzError , it->second.second.dxdz() , it->second.second.dxdzError());
368  weight(dydz , dydzError , it->second.second.dydz() , it->second.second.dydzError());
369  weight(widthX, widthXError, it->second.second.BeamWidthX(), it->second.second.BeamWidthXError());
370  weight(widthY, widthYError, it->second.second.BeamWidthY(), it->second.second.BeamWidthYError());
371  if(it->second.second.type() == BeamSpot::Tracker){
372  type = BeamSpot::Tracker;
373  }
374  }
375  BeamSpot::Point bsPosition(x,y,z);
377  error(0,0) = xError*xError;
378  error(1,1) = yError*yError;
379  error(2,2) = zError*zError;
380  error(3,3) = sigmaZError*sigmaZError;
381  error(4,4) = dxdzError*dxdzError;
382  error(5,5) = dydzError*dydzError;
383  error(6,6) = widthXError*widthXError;
384  BeamSpot weightedBeamSpot(bsPosition,sigmaZ,dxdz,dydz,widthX,error,type);
385  weightedBeamSpot.setBeamWidthY(widthY);
386  LogInfo("AlcaBeamSpotManager")
387  << "Weighted BeamSpot will be:" <<'\n'
388  << weightedBeamSpot
389  << endl;
390  return weightedBeamSpot;
391 }
392 
393 //--------------------------------------------------------------------------------------------------
394 void AlcaBeamSpotManager::weight(double& mean,double& meanError,const double& val,const double& valError){
395  double tmpError = 0;
396  if (meanError < 1e-8){
397  tmpError = 1/(valError*valError);
398  mean = val*tmpError;
399  }
400  else{
401  tmpError = 1/(meanError*meanError) + 1/(valError*valError);
402  mean = mean/(meanError*meanError) + val/(valError*valError);
403  }
404  mean = mean/tmpError;
405  meanError = sqrt(1/tmpError);
406 }
407 
408 //--------------------------------------------------------------------------------------------------
409 pair<float,float> AlcaBeamSpotManager::delta(const float& x, const float& xError, const float& nextX, const float& nextXError){
410  return pair<float,float>(x - nextX, sqrt(pow(xError,2) + pow(nextXError,2)) );
411 }
412 
413 //--------------------------------------------------------------------------------------------------
414 float AlcaBeamSpotManager::deltaSig(const float& num, const float& den){
415  if(den != 0){
416  return fabs(num/den);
417  }
418  else{
419  return float(LONG_MAX);
420  }
421 }
422 
type
Definition: HCALResponse.h:21
math::Error< dimension >::type CovarianceMatrix
Definition: BeamSpot.h:31
double z0() const
z coordinate
Definition: BeamSpot.h:68
BeamType
beam spot flags
Definition: BeamSpot.h:26
double sigmaZ0Error() const
error on sigma z
Definition: BeamSpot.h:96
std::string beamSpotOutputBase_
edm::InputTag beamSpotTag_
edm::EDGetTokenT< reco::BeamSpot > beamSpotToken_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
double dydzError() const
error on dydz
Definition: BeamSpot.h:100
Definition: weight.py:1
math::XYZPoint Point
point in the space
Definition: BeamSpot.h:29
Timestamp const & beginTime() const
std::map< edm::LuminosityBlockNumber_t, std::pair< edm::Timestamp, reco::BeamSpot > > beamSpotMap_
unsigned int LuminosityBlockNumber_t
std::map< edm::LuminosityBlockNumber_t, std::pair< edm::Timestamp, reco::BeamSpot > >::iterator bsMap_iterator
std::string beamSpotModuleName_
LuminosityBlockNumber_t luminosityBlock() const
double dydz() const
dydz slope
Definition: BeamSpot.h:84
void setBeamWidthY(double v)
Definition: BeamSpot.h:109
double dxdzError() const
error on dxdz
Definition: BeamSpot.h:98
T sqrt(T t)
Definition: SSEVec.h:18
reco::BeamSpot weight(const bsMap_iterator &begin, const bsMap_iterator &end)
double BeamWidthX() const
beam width X
Definition: BeamSpot.h:86
virtual ~AlcaBeamSpotManager(void)
#define end
Definition: vmac.h:39
double BeamWidthYError() const
error on beam width Y, assume error in X = Y
Definition: BeamSpot.h:105
bool isValid() const
Definition: HandleBase.h:74
double BeamWidthXError() const
error on beam width X, assume error in X = Y
Definition: BeamSpot.h:103
double z0Error() const
error on z
Definition: BeamSpot.h:94
double dxdz() const
dxdz slope
Definition: BeamSpot.h:82
double x0Error() const
error on x
Definition: BeamSpot.h:90
double y0Error() const
error on y
Definition: BeamSpot.h:92
void readLumi(const edm::LuminosityBlock &)
T const * product() const
Definition: Handle.h:81
double sigmaZ() const
sigma z
Definition: BeamSpot.h:80
double BeamWidthY() const
beam width Y
Definition: BeamSpot.h:88
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
double y0() const
y coordinate
Definition: BeamSpot.h:66
std::pair< float, float > delta(const float &x, const float &xError, const float &nextX, const float &nextXError)
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
float deltaSig(const float &num, const float &den)
double x0() const
x coordinate
Definition: BeamSpot.h:64