CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TkLasBeamFitter.cc
Go to the documentation of this file.
1 
17 // framework include files
29 
35 
36 // data formats
37 // for edm::InRun
39 
40 // laser data formats
43 
44 // further includes
49 
50 #include <iostream>
51 #include "TMinuit.h"
52 #include "TGraphErrors.h"
53 #include "TF1.h"
54 #include "TH1.h"
55 #include "TH2.h"
56 
57 using namespace edm;
58 using namespace std;
59 
60 
61 //
62 // class declaration
63 //
64 
65 class TkLasBeamFitter : public edm::one::EDProducer<edm::EndRunProducer> {
66 public:
67  explicit TkLasBeamFitter(const edm::ParameterSet &config);
68  ~TkLasBeamFitter();
69 
70  //virtual void beginJob(const edm::EventSetup& /*access deprecated*/) {}
71  virtual void produce(edm::Event &event, const edm::EventSetup &setup) override;
72  // virtual void beginRun(edm::Run &run, const edm::EventSetup &setup);
73  virtual void endRunProduce(edm::Run &run, const edm::EventSetup &setup) override;
74  //virtual void endJob() {}
75 
76 private:
79  void getLasBeams(TkFittedLasBeam &beam,vector<TrajectoryStateOnSurface> &tsosLas);
80  void getLasHits(TkFittedLasBeam &beam, const SiStripLaserRecHit2D &hit,
81  vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globHit,
82  unsigned int &hitsAtTecPlus);
83 // void fillVectors(TkFittedLasBeam &beam);
84 
85  // need static functions to be used in fitter(..);
86  // all parameters used therein have to be static, as well (see below)
87  static double tecPlusFunction(double *x, double *par);
88  static double tecMinusFunction(double *x, double *par);
89  static double atFunction(double *x, double *par);
90 
91  void fitter(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix,
92  unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
93  std::vector<double> &hitPhi, std::vector<double> &hitPhiError, std::vector<double> &hitZprimeError,
94  double &zMin, double &zMax, double &bsAngleParam,
95  double &offset, double &offsetError, double &slope, double &slopeError,
96  double &phiAtMinusParam, double &phiAtPlusParam,
97  double &atThetaSplitParam);
98 
99  void trackPhi(TkFittedLasBeam &beam, unsigned int &hit,
100  double &trackPhi, double &trackPhiRef,
101  double &offset, double &slope, double &bsAngleParam,
102  double &phiAtMinusParam, double &phiAtPlusParam,
103  double &atThetaSplitParam, std::vector<GlobalPoint> &globHit);
104 
105  void globalTrackPoint(TkFittedLasBeam &beam,
106  unsigned int &hit, unsigned int &hitsAtTecPlus,
107  double &trackPhi, double &trackPhiRef,
108  std::vector<GlobalPoint> &globHit, std::vector<GlobalPoint> &globPtrack,
109  GlobalPoint &globPref, std::vector<double> &hitPhiError);
110 
111  void buildTrajectory(TkFittedLasBeam &beam, unsigned int &hit,
112  vector<const GeomDetUnit*> &gd, std::vector<GlobalPoint> &globPtrack,
113  vector<TrajectoryStateOnSurface> &tsosLas, GlobalPoint &globPref);
114 
115  bool fitBeam(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix,
116  unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
117  double &offset, double &slope, vector<GlobalPoint> &globPtrack,
118  double &bsAngleParam, double &chi2);
119 
120  // ----------member data ---------------------------
123  unsigned int nAtParameters_;
124 
126 
127  // static parameters used in static parametrization functions
128  static vector<double> gHitZprime;
129  static vector<double> gBarrelModuleRadius;
130  static vector<double> gBarrelModuleOffset;
131  static float gTIBparam;
132  static float gTOBparam;
133  static double gBeamR;
134  static double gBeamZ0;
135  static double gBeamSplitterZprime;
136  static unsigned int gHitsAtTecMinus;
137  static double gBSparam;
138  static bool gFitBeamSplitters;
139  static bool gIsInnerBarrel;
140 
141  // histograms
142  TH1F *h_bsAngle, *h_hitX, *h_hitXTecPlus, *h_hitXTecMinus,
143  *h_hitXAt, *h_chi2, *h_chi2ndof, *h_pull, *h_res,
144  *h_resTecPlus, *h_resTecMinus, *h_resAt;
145  TH2F *h_bsAngleVsBeam, *h_hitXvsZTecPlus, *h_hitXvsZTecMinus,
146  *h_hitXvsZAt, *h_resVsZTecPlus, *h_resVsZTecMinus, *h_resVsZAt,
147  *h_resVsHitTecPlus, *h_resVsHitTecMinus, *h_resVsHitAt;
148 };
149 
150 //
151 // constants, enums and typedefs
152 //
153 
154 
155 //
156 // static data member definitions
157 //
158 
159 // static parameters used in parametrization functions
160 vector<double> TkLasBeamFitter::gHitZprime;
163 float TkLasBeamFitter::gTIBparam = 0.097614; // = abs(r_offset/r_module) (nominal!)
164 float TkLasBeamFitter::gTOBparam = 0.034949; // = abs(r_offset/r_module) (nominal!)
165 double TkLasBeamFitter::gBeamR = 0.0;
166 double TkLasBeamFitter::gBeamZ0 = 0.0;
168 unsigned int TkLasBeamFitter::gHitsAtTecMinus = 0;
169 double TkLasBeamFitter::gBSparam = 0.0;
172 
173 // handles
177 
178 //
179 // constructors and destructor
180 //
182  src_(iConfig.getParameter<edm::InputTag>("src")),
183  fitBeamSplitters_(iConfig.getParameter<bool>("fitBeamSplitters")),
184  nAtParameters_(iConfig.getParameter<unsigned int>("numberOfFittedAtParameters")),
185  h_bsAngle(0), h_hitX(0), h_hitXTecPlus(0), h_hitXTecMinus(0),
186  h_hitXAt(0), h_chi2(0), h_chi2ndof(0), h_pull(0), h_res(0),
187  h_resTecPlus(0), h_resTecMinus(0), h_resAt(0),
188  h_bsAngleVsBeam(0), h_hitXvsZTecPlus(0), h_hitXvsZTecMinus(0),
189  h_hitXvsZAt(0), h_resVsZTecPlus(0), h_resVsZTecMinus(0), h_resVsZAt(0),
190  h_resVsHitTecPlus(0), h_resVsHitTecMinus(0), h_resVsHitAt(0)
191 {
192  // declare the products to produce
193  this->produces<TkFittedLasBeamCollection, edm::InRun>();
194  this->produces<TsosVectorCollection, edm::InRun>();
195 
196  //now do what ever other initialization is needed
197 }
198 
199 //---------------------------------------------------------------------------------------
201 {
202 
203  // do anything here that needs to be done at desctruction time
204  // (e.g. close files, deallocate resources etc.)
205 
206 }
207 
208 
209 //
210 // member functions
211 //
212 
213 //---------------------------------------------------------------------------------------
214 // ------------ method called to produce the data ------------
216 {
217  // Nothing per event!
218 }
219 
220 //---------------------------------------------------------------------------------------
221 // ------------ method called at end of each run ---------------------------------------
223 {
224 // }
225 // // FIXME!
226 // // Indeed, that should be in endRun(..) - as soon as AlignmentProducer can call
227 // // the algorithm's endRun correctly!
228 //
229 //
230 // void TkLasBeamFitter::beginRun(edm::Run &run, const edm::EventSetup &setup)
231 // {
232 
233  // book histograms
234  h_hitX = fs->make<TH1F>("hitX","local x of LAS hits;local x [cm];N",100,-0.5,0.5);
235  h_hitXTecPlus = fs->make<TH1F>("hitXTecPlus","local x of LAS hits in TECplus;local x [cm];N",100,-0.5,0.5);
236  h_hitXTecMinus = fs->make<TH1F>("hitXTecMinus","local x of LAS hits in TECminus;local x [cm];N",100,-0.5,0.5);
237  h_hitXAt = fs->make<TH1F>("hitXAt","local x of LAS hits in ATs;local x [cm];N",100,-2.5,2.5);
238  h_hitXvsZTecPlus = fs->make<TH2F>("hitXvsZTecPlus","local x vs z in TECplus;z [cm];local x [cm]",80,120,280,100,-0.5,0.5);
239  h_hitXvsZTecMinus = fs->make<TH2F>("hitXvsZTecMinus","local x vs z in TECMinus;z [cm];local x [cm]",80,-280,-120,100,-0.5,0.5);
240  h_hitXvsZAt = fs->make<TH2F>("hitXvsZAt","local x vs z in ATs;z [cm];local x [cm]",200,-200,200,100,-0.5,0.5);
241  h_chi2 = fs->make<TH1F>("chi2","#chi^{2};#chi^{2};N",100,0,2000);
242  h_chi2ndof = fs->make<TH1F>("chi2ndof","#chi^{2} per degree of freedom;#chi^{2}/N_{dof};N",100,0,300);
243  h_pull = fs->make<TH1F>("pull","pulls of #phi residuals;pull;N",50,-10,10);
244  h_res = fs->make<TH1F>("res","#phi residuals;#phi_{track} - #phi_{hit} [rad];N",60,-0.0015,0.0015);
245  h_resTecPlus = fs->make<TH1F>("resTecPlus","#phi residuals in TECplus;#phi_{track} - #phi_{hit} [rad];N",30,-0.0015,0.0015);
246  h_resTecMinus = fs->make<TH1F>("resTecMinus","#phi residuals in TECminus;#phi_{track} - #phi_{hit} [rad];N",60,-0.0015,0.0015);
247  h_resAt = fs->make<TH1F>("resAt","#phi residuals in ATs;#phi_{track} - #phi_{hit} [rad];N",30,-0.0015,0.0015);
248  h_resVsZTecPlus = fs->make<TH2F>("resVsZTecPlus","phi residuals vs. z in TECplus;z [cm];#phi_{track} - #phi_{hit} [rad]",
249  80,120,280,100,-0.0015,0.0015);
250  h_resVsZTecMinus = fs->make<TH2F>("resVsZTecMinus","phi residuals vs. z in TECminus;z [cm];#phi_{track} - #phi_{hit} [rad]",
251  80,-280,-120,100,-0.0015,0.0015);
252  h_resVsZAt = fs->make<TH2F>("resVsZAt","#phi residuals vs. z in ATs;N;#phi_{track} - #phi_{hit} [rad]",
253  200,-200,200,100,-0.0015,0.0015);
254  h_resVsHitTecPlus = fs->make<TH2F>("resVsHitTecPlus","#phi residuals vs. hits in TECplus;hit no.;#phi_{track} - #phi_{hit} [rad]",
255  144,0,144,100,-0.0015,0.0015);
256  h_resVsHitTecMinus = fs->make<TH2F>("resVsHitTecMinus","#phi residuals vs. hits in TECminus;hit no.;#phi_{track} - #phi_{hit} [rad]",
257  144,0,144,100,-0.0015,0.0015);
258  h_resVsHitAt = fs->make<TH2F>("resVsHitAt","#phi residuals vs. hits in ATs;hit no.;#phi_{track} - #phi_{hit} [rad]",
259  176,0,176,100,-0.0015,0.0015);
260  h_bsAngle = fs->make<TH1F>("bsAngle","fitted beam splitter angle;BS angle [rad];N",40,-0.004,0.004);
261  h_bsAngleVsBeam = fs->make<TH2F>("bsAngleVsBeam","fitted beam splitter angle per beam;Beam no.;BS angle [rad]",
262  40,0,300,100,-0.004,0.004);
263 
264  // Create output collections - they are parallel.
265  // (edm::Ref etc. and thus edm::AssociationVector are not supported for edm::Run...)
266  std::auto_ptr<TkFittedLasBeamCollection> fittedBeams(new TkFittedLasBeamCollection);
267  // One std::vector<TSOS> for each TkFittedLasBeam:
268  std::auto_ptr<TsosVectorCollection> tsosesVec(new TsosVectorCollection);
269 
270  // get TkLasBeams, Tracker geometry, magnetic field
271  run.getByLabel( "LaserAlignment", "tkLaserBeams", laserBeams );
272  setup.get<TrackerDigiGeometryRecord>().get(geometry);
274 
275  // hack for fixed BSparams (ugly!)
276 // double bsParams[34] = {-0.000266,-0.000956,-0.001205,-0.000018,-0.000759,0.002554,
277 // 0.000465,0.000975,0.001006,0.002027,-0.001263,-0.000763,
278 // -0.001702,0.000906,-0.002120,0.001594,0.000661,-0.000457,
279 // -0.000447,0.000347,-0.002266,-0.000446,0.000659,0.000018,
280 // -0.001630,-0.000324,
281 // // ATs
282 // -999.,-0.001709,-0.002091,-999.,
283 // -0.001640,-999.,-0.002444,-0.002345};
284 
285  double bsParams[40] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,
286  0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0};
287 
288  // beam counter
289  unsigned int beamNo(0);
290  // fit BS? If false, values from bsParams are taken
292  if(fitBeamSplitters_) cout << "Fitting BS!" << endl;
293  else cout << "BS fixed, not fitted!" << endl;
294 
295  // loop over LAS beams
296  for(TkLasBeamCollection::const_iterator iBeam = laserBeams->begin(), iEnd = laserBeams->end();
297  iBeam != iEnd; ++iBeam){
298 
299  TkFittedLasBeam beam(*iBeam);
300  vector<TrajectoryStateOnSurface> tsosLas;
301 
302  // set BS param for fit
303  gBSparam = bsParams[beamNo];
304 
305  // call main function; all other functions are called inside getLasBeams(..)
306  this->getLasBeams(beam, tsosLas);
307 
308  // fill output products
309  fittedBeams->push_back(beam);
310  tsosesVec->push_back(tsosLas);
311 
312 // if(!this->fitBeam(fittedBeams->back(), tsosesVec->back())){
313 // edm::LogError("BadFit")
314 // << "Problems fitting TkLasBeam, id " << fittedBeams->back().getBeamId() << ".";
315 // fittedBeams->pop_back(); // remove last entry added just before
316 // tsosesVec->pop_back(); // dito
317 // }
318 
319  beamNo++;
320  }
321 
322  // finally, put fitted beams and TSOS vectors into run
323  run.put(fittedBeams);
324  run.put(tsosesVec);
325 }
326 
327 // methods for las data processing
328 
329 // -------------- loop over beams, call functions ----------------------------
330 void TkLasBeamFitter::getLasBeams(TkFittedLasBeam &beam, vector<TrajectoryStateOnSurface> &tsosLas)
331 {
332  cout << "---------------------------------------" << endl;
333  cout << "beam id: " << beam.getBeamId() // << " isTec: " << (beam.isTecInternal() ? "Y" : "N")
334  << " isTec+: " << (beam.isTecInternal(1) ? "Y" : "N") << " isTec-: " << (beam.isTecInternal(-1) ? "Y" : "N")
335  << " isAt: " << (beam.isAlignmentTube() ? "Y" : "N") << " isR6: " << (beam.isRing6() ? "Y" : "N")
336  << endl;
337 
338  // reset static variables
339  gHitsAtTecMinus = 0;
340  gHitZprime.clear();
341  gBarrelModuleRadius.clear();
342  gBarrelModuleOffset.clear();
343 
344  // set right beam radius
345  gBeamR = beam.isRing6() ? 84.0 : 56.4;
346 
347  vector<const GeomDetUnit*> gd;
348  vector<GlobalPoint> globHit;
349  unsigned int hitsAtTecPlus(0);
350  double sumZ(0.);
351 
352  // loop over hits
353  for( TkLasBeam::const_iterator iHit = beam.begin(); iHit < beam.end(); ++iHit ){
354  // iHit is a SiStripLaserRecHit2D
355 
356  const SiStripLaserRecHit2D hit(*iHit);
357 
358  this->getLasHits(beam, hit, gd, globHit, hitsAtTecPlus);
359  sumZ += globHit.back().z();
360 
361  // fill histos
362  h_hitX->Fill(hit.localPosition().x());
363  // TECplus
364  if(beam.isTecInternal(1)){
365  h_hitXTecPlus->Fill(hit.localPosition().x());
366  h_hitXvsZTecPlus->Fill(globHit.back().z(),hit.localPosition().x());
367  }
368  // TECminus
369  else if(beam.isTecInternal(-1)){
370  h_hitXTecMinus->Fill(hit.localPosition().x());
371  h_hitXvsZTecMinus->Fill(globHit.back().z(),hit.localPosition().x());
372  }
373  // ATs
374  else{
375  h_hitXAt->Fill(hit.localPosition().x());
376  h_hitXvsZAt->Fill(globHit.back().z(),hit.localPosition().x());
377  }
378  }
379 
380  gBeamZ0 = sumZ / globHit.size();
381  double zMin(0.), zMax(0.);
382  // TECplus
383  if(beam.isTecInternal(1)){
384  gBeamSplitterZprime = 205.75 - gBeamZ0;
385  zMin = 120.0 - gBeamZ0;
386  zMax = 280.0 - gBeamZ0;
387  }
388  // TECminus
389  else if(beam.isTecInternal(-1)){
390  gBeamSplitterZprime = -205.75 - gBeamZ0;
391  zMin = -280.0 - gBeamZ0;
392  zMax = -120.0 - gBeamZ0;
393  }
394  // AT
395  else{
396  gBeamSplitterZprime = 112.3 - gBeamZ0;
397  zMin = -200.0 - gBeamZ0;
398  zMax = 200.0 - gBeamZ0;
399  }
400 
401  // fill vectors for fitted quantities
402  vector<double> hitPhi, hitPhiError, hitZprimeError;
403 
404  for(unsigned int hit = 0; hit < globHit.size(); ++hit){
405  hitPhi.push_back(static_cast<double>(globHit[hit].phi()));
406  // localPositionError[hit] or assume 0.003, 0.006
407  hitPhiError.push_back( 0.003 / globHit[hit].perp());
408  // no errors on z, fill with zeros
409  hitZprimeError.push_back(0.0);
410  // barrel-specific values
411  if(beam.isAlignmentTube() && abs(globHit[hit].z()) < 112.3){
412  gBarrelModuleRadius.push_back(globHit[hit].perp());
413  gBarrelModuleOffset.push_back(gBarrelModuleRadius.back() - gBeamR);
414  // TIB/TOB flag
415  if(gBarrelModuleOffset.back() < 0.0){
416  gIsInnerBarrel = 1;
417  }
418  else{
419  gIsInnerBarrel = 0;
420  }
421  gHitZprime.push_back(globHit[hit].z() - gBeamZ0 - abs(gBarrelModuleOffset.back()));
422  }
423  // non-barrel z'
424  else{
425  gHitZprime.push_back(globHit[hit].z() - gBeamZ0);
426  }
427  }
428 
429  // number of fit parameters, 3 for TECs (always!); 3, 5, or 6 for ATs
430  unsigned int tecParams(3), atParams(0);
431  if(nAtParameters_ == 3) atParams = 3;
432  else if(nAtParameters_ == 5) atParams = 5;
433  else atParams = 6; // <-- default value, recommended
434  unsigned int nFitParams(0);
435  if(!fitBeamSplitters_ ||
436  (hitsAtTecPlus == 0 && beam.isAlignmentTube() ) ){
437  tecParams = tecParams - 1;
438  atParams = atParams - 1;
439  }
440  if(beam.isTecInternal()){
441  nFitParams = tecParams;
442  }
443  else{
444  nFitParams = atParams;
445  }
446 
447  // fit parameter definitions
448  double offset(0.), offsetError(0.), slope(0.), slopeError(0.),
449  bsAngleParam(0.), phiAtMinusParam(0.), phiAtPlusParam(0.),
450  atThetaSplitParam(0.);
451  AlgebraicSymMatrix covMatrix;
452  if(!fitBeamSplitters_ || (beam.isAlignmentTube() && hitsAtTecPlus == 0)){
453  covMatrix = AlgebraicSymMatrix(nFitParams, 1);
454  }
455  else{
456  covMatrix = AlgebraicSymMatrix(nFitParams - 1, 1);
457  }
458 
459  this->fitter(beam, covMatrix,
460  hitsAtTecPlus, nFitParams,
461  hitPhi, hitPhiError, hitZprimeError,
462  zMin, zMax, bsAngleParam,
463  offset, offsetError, slope, slopeError,
464  phiAtMinusParam, phiAtPlusParam,
465  atThetaSplitParam);
466 
467  vector<GlobalPoint> globPtrack;
468  GlobalPoint globPref;
469  double chi2(0.);
470 
471  for(unsigned int hit = 0; hit < gHitZprime.size(); ++hit){
472 
473  // additional phi value (trackPhiRef) for trajectory calculation
474  double trackPhi(0.), trackPhiRef(0.);
475 
476  this->trackPhi(beam, hit, trackPhi, trackPhiRef,
477  offset, slope, bsAngleParam,
478  phiAtMinusParam, phiAtPlusParam,
479  atThetaSplitParam, globHit);
480 
481  cout << "track phi = " << trackPhi
482  << ", hit phi = " << hitPhi[hit]
483  << ", zPrime = " << gHitZprime[hit]
484  << " r = " << globHit[hit].perp() << endl;
485 
486  this->globalTrackPoint(beam, hit, hitsAtTecPlus,
487  trackPhi, trackPhiRef,
488  globHit, globPtrack, globPref,
489  hitPhiError);
490 
491  // calculate residuals = pred - hit (in global phi)
492  const double phiResidual = globPtrack[hit].phi() - globHit[hit].phi();
493  // pull calculation (FIX!!!)
494  const double phiResidualPull = phiResidual / hitPhiError[hit];
495  // sqrt(hitPhiError[hit]*hitPhiError[hit] +
496  // (offsetError*offsetError + globPtrack[hit].z()*globPtrack[hit].z() * slopeError*slopeError));
497 
498  // calculate chi2
499  chi2 += phiResidual*phiResidual / (hitPhiError[hit]*hitPhiError[hit]);
500 
501  // fill histos
502  h_res->Fill(phiResidual);
503  // TECplus
504  if(beam.isTecInternal(1)){
505  h_pull->Fill(phiResidualPull);
506  h_resTecPlus->Fill(phiResidual);
507  h_resVsZTecPlus->Fill(globPtrack[hit].z(), phiResidual);
508  // Ring 6
509  if(beam.isRing6()){
510  h_resVsHitTecPlus->Fill(hit+(beam.getBeamId()-1)/10*9+72, phiResidual);
511  }
512  // Ring 4
513  else{
514  h_resVsHitTecPlus->Fill(hit+beam.getBeamId()/10*9, phiResidual);
515  }
516  }
517  // TECminus
518  else if(beam.isTecInternal(-1)){
519  h_pull->Fill(phiResidualPull);
520  h_resTecMinus->Fill(phiResidual);
521  h_resVsZTecMinus->Fill(globPtrack[hit].z(), phiResidual);
522  // Ring 6
523  if(beam.isRing6()){
524  h_resVsHitTecMinus->Fill(hit+(beam.getBeamId()-101)/10*9+72, phiResidual);
525  }
526  // Ring 4
527  else{
528  h_resVsHitTecMinus->Fill(hit+(beam.getBeamId()-100)/10*9, phiResidual);
529  }
530  }
531  // ATs
532  else{
533  h_pull->Fill(phiResidualPull);
534  h_resAt->Fill(phiResidual);
535  h_resVsZAt->Fill(globPtrack[hit].z(), phiResidual);
536  h_resVsHitAt->Fill(hit+(beam.getBeamId()-200)/10*22, phiResidual);
537  }
538 
539  this->buildTrajectory(beam, hit, gd, globPtrack, tsosLas, globPref);
540  }
541 
542  cout << "chi^2 = " << chi2 << ", chi^2/ndof = " << chi2/(gHitZprime.size() - nFitParams) << endl;
543  this->fitBeam(beam, covMatrix, hitsAtTecPlus, nFitParams,
544  offset, slope, globPtrack, bsAngleParam, chi2);
545 
546  cout << "bsAngleParam = " << bsAngleParam << endl;
547 
548  // fill histos
549  // include slope, offset, covariance plots here
550  h_chi2->Fill(chi2);
551  h_chi2ndof->Fill(chi2/(gHitZprime.size() - nFitParams));
552  if(bsAngleParam != 0.0){
553  h_bsAngle->Fill(2.0*atan(0.5*bsAngleParam));
554  h_bsAngleVsBeam->Fill(beam.getBeamId(), 2.0*atan(0.5*bsAngleParam));
555  }
556 }
557 
558 
559 // --------- get hits, convert to global coordinates ---------------------------
561  vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globHit,
562  unsigned int &hitsAtTecPlus)
563 {
564  // get global position of LAS hits
565  gd.push_back(geometry->idToDetUnit(hit.getDetId()));
566  GlobalPoint globPtemp(gd.back()->toGlobal(hit.localPosition()));
567 
568  // testing: globPtemp should be right
569  globHit.push_back(globPtemp);
570 
571  if(beam.isAlignmentTube()){
572  if(abs(globPtemp.z()) > 112.3){
573  if(globPtemp.z() < 112.3) gHitsAtTecMinus++ ;
574  else hitsAtTecPlus++ ;
575  }
576  }
577 }
578 
579 
580 // ------------ parametrization functions for las beam fits ------------
581 double TkLasBeamFitter::tecPlusFunction(double *x, double *par)
582 {
583  double z = x[0]; // 'primed'? -> yes!!!
584 
585  if(z < gBeamSplitterZprime){
586  return par[0] + par[1] * z;
587  }
588  else{
589  if(gFitBeamSplitters){
590  // par[2] = 2*tan(BeamSplitterAngle/2.0)
591  return par[0] + par[1] * z - par[2] * (z - gBeamSplitterZprime)/gBeamR;
592  }
593  else{
594  return par[0] + par[1] * z - gBSparam * (z - gBeamSplitterZprime)/gBeamR;
595  }
596  }
597 }
598 
599 
600 double TkLasBeamFitter::tecMinusFunction(double *x, double *par)
601 {
602  double z = x[0]; // 'primed'? -> yes!!!
603 
604  if(z > gBeamSplitterZprime){
605  return par[0] + par[1] * z;
606  }
607  else{
608  if(gFitBeamSplitters){
609  // par[2] = 2*tan(BeamSplitterAngle/2.0)
610  return par[0] + par[1] * z + par[2] * (z - gBeamSplitterZprime)/gBeamR;
611  }
612  else{
613  return par[0] + par[1] * z + gBSparam * (z - gBeamSplitterZprime)/gBeamR;
614  }
615  }
616 }
617 
618 double TkLasBeamFitter::atFunction(double *x, double *par)
619 {
620  double z = x[0]; // 'primed'? -> yes!!!
621  // TECminus
622  if(z < -gBeamSplitterZprime - 2.0*gBeamZ0){
623  return par[0] + par[1] * z;
624  }
625  // BarrelMinus
626  else if(-gBeamSplitterZprime - 2.0*gBeamZ0 < z && z < -gBeamZ0){
627  // z value includes module offset from main beam axis
628  // TOB
629  if(!gIsInnerBarrel){
630  return par[0] + par[1] * z + gTOBparam * (par[2] + par[4]);
631  }
632  // TIB
633  else{
634  return par[0] + par[1] * z - gTIBparam * (par[2] - par[4]);
635  }
636  }
637  // BarrelPlus
638  else if(-gBeamZ0 < z && z < gBeamSplitterZprime){
639  // z value includes module offset from main beam axis
640  // TOB
641  if(!gIsInnerBarrel){
642  return par[0] + par[1] * z + gTOBparam * (par[3] - par[4]);
643  }
644  // TIB
645  else{
646  return par[0] + par[1] * z - gTIBparam * (par[3] + par[4]);
647  }
648  }
649  // TECplus
650  else{
651  if(gFitBeamSplitters){
652  // par[2] = 2*tan(BeamSplitterAngle/2.0)
653  return par[0] + par[1] * z - par[5] * (z - gBeamSplitterZprime)/gBeamR; // BS par: 5, 4, or 2
654  }
655  else{
656  return par[0] + par[1] * z - gBSparam * (z - gBeamSplitterZprime)/gBeamR;
657  }
658  }
659 }
660 
661 
662 // ------------ perform fit of beams ------------------------------------
664  unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
665  vector<double> &hitPhi, vector<double> &hitPhiError, vector<double> &hitZprimeError,
666  double &zMin, double &zMax, double &bsAngleParam,
667  double &offset, double &offsetError, double &slope, double &slopeError,
668  double &phiAtMinusParam, double &phiAtPlusParam, double &atThetaSplitParam)
669 {
670  TGraphErrors *lasData = new TGraphErrors(gHitZprime.size(),
671  &(gHitZprime[0]), &(hitPhi[0]),
672  &(hitZprimeError[0]), &(hitPhiError[0]));
673 
674  // do fit (R = entire range)
675  if(beam.isTecInternal(1)){
676  TF1 tecPlus("tecPlus", tecPlusFunction, zMin, zMax, nFitParams );
677  tecPlus.SetParameter( 1, 0 ); // slope
678  tecPlus.SetParameter( nFitParams - 1, 0 ); // BS
679  lasData->Fit(&tecPlus, "R"); // "R", "RV" or "RQ"
680  }
681  else if(beam.isTecInternal(-1)){
682  TF1 tecMinus("tecMinus", tecMinusFunction, zMin, zMax, nFitParams );
683  tecMinus.SetParameter( 1, 0 ); // slope
684  tecMinus.SetParameter( nFitParams - 1, 0 ); // BS
685  lasData->Fit(&tecMinus, "R");
686  }
687  else{
688  TF1 at("at", atFunction, zMin, zMax, nFitParams );
689  at.SetParameter( 1, 0 ); // slope
690  at.SetParameter( nFitParams - 1, 0 ); // BS
691  lasData->Fit(&at,"R");
692  }
693 
694  // get values and errors for offset and slope
695  gMinuit->GetParameter(0, offset, offsetError);
696  gMinuit->GetParameter(1, slope, slopeError);
697 
698  // additional AT parameters
699  // define param errors that are not used later
700  double bsAngleParamError(0.), phiAtMinusParamError(0.),
701  phiAtPlusParamError(0.), atThetaSplitParamError(0.);
702 
703  if(beam.isAlignmentTube()){
704  gMinuit->GetParameter(2, phiAtMinusParam, phiAtMinusParamError);
705  gMinuit->GetParameter(3, phiAtPlusParam, phiAtPlusParamError);
706  gMinuit->GetParameter(4, atThetaSplitParam, atThetaSplitParamError);
707  }
708  // get Beam Splitter parameters
709  if(fitBeamSplitters_){
710  if(beam.isAlignmentTube() && hitsAtTecPlus == 0){
711  bsAngleParam = gBSparam;
712  }
713  else{
714  gMinuit->GetParameter( nFitParams - 1 , bsAngleParam, bsAngleParamError);
715  }
716  }
717  else{
718  bsAngleParam = gBSparam;
719  }
720 
721  // fill covariance matrix
722  vector<double> vec( covMatrix.num_col() * covMatrix.num_col() );
723  gMinuit->mnemat( &vec[0], covMatrix.num_col() );
724  for(int col = 0; col < covMatrix.num_col(); col++){
725  for(int row = 0; row < covMatrix.num_col(); row++){
726  covMatrix[col][row] = vec[row + covMatrix.num_col()*col];
727  }
728  }
729  // compute correlation between parameters
730 // double corr01 = covMatrix[1][0]/(offsetError*slopeError);
731 
732  delete lasData;
733 }
734 
735 
736 // -------------- calculate track phi value ----------------------------------
738  double &trackPhi, double &trackPhiRef,
739  double &offset, double &slope, double &bsAngleParam,
740  double &phiAtMinusParam, double &phiAtPlusParam,
741  double &atThetaSplitParam, vector<GlobalPoint> &globHit)
742 {
743  // TECplus
744  if(beam.isTecInternal(1)){
745  if(gHitZprime[hit] < gBeamSplitterZprime){
746  trackPhi = offset + slope * gHitZprime[hit];
747  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
748  }
749  else{
750  trackPhi = offset + slope * gHitZprime[hit]
751  - bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
752  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
753  - bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
754  }
755  }
756  // TECminus
757  else if(beam.isTecInternal(-1)){
758  if(gHitZprime[hit] > gBeamSplitterZprime){
759  trackPhi = offset + slope * gHitZprime[hit];
760  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
761  }
762  else{
763  trackPhi = offset + slope * gHitZprime[hit]
764  + bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
765  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
766  + bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
767  }
768  }
769  // ATs
770  else{
771  // TECminus
772  if(gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0){
773  trackPhi = offset + slope * gHitZprime[hit];
774  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0);
775  }
776  // BarrelMinus
777  else if(-gBeamSplitterZprime - 2.0*gBeamZ0 < gHitZprime[hit] && gHitZprime[hit] < -gBeamZ0){
778  if(!gIsInnerBarrel){
779  trackPhi = offset + slope * gHitZprime[hit] + gTOBparam * (phiAtMinusParam + atThetaSplitParam);
780  }
781  else{
782  trackPhi = offset + slope * gHitZprime[hit] - gTIBparam * (phiAtMinusParam - atThetaSplitParam);
783  }
784  trackPhiRef = offset + slope * (gHitZprime[hit] + abs(gBarrelModuleOffset[hit-gHitsAtTecMinus]));
785  }
786  // BarrelPlus
787  else if(-gBeamZ0 < gHitZprime[hit] && gHitZprime[hit] < gBeamSplitterZprime){
788  if(!gIsInnerBarrel){
789  trackPhi = offset + slope * gHitZprime[hit] + gTOBparam * (phiAtPlusParam - atThetaSplitParam);
790  }
791  else{
792  trackPhi = offset + slope * gHitZprime[hit] - gTIBparam * (phiAtPlusParam + atThetaSplitParam);
793  }
794  trackPhiRef = offset + slope * (gHitZprime[hit] + abs(gBarrelModuleOffset[hit-gHitsAtTecMinus]));
795  }
796  // TECplus
797  else{
798  trackPhi = offset + slope * gHitZprime[hit]
799  - bsAngleParam * (gHitZprime[hit] - gBeamSplitterZprime)/gBeamR;
800  trackPhiRef = offset + slope * (gHitZprime[hit] + 1.0)
801  - bsAngleParam * ((gHitZprime[hit] + 1.0) - gBeamSplitterZprime)/gBeamR;
802  }
803  }
804 }
805 
806 
807 // -------------- calculate global track points, hit residuals, chi2 ----------------------------------
809  unsigned int &hit, unsigned int &hitsAtTecPlus,
810  double &trackPhi, double &trackPhiRef,
811  vector<GlobalPoint> &globHit, vector<GlobalPoint> &globPtrack,
812  GlobalPoint &globPref, vector<double> &hitPhiError)
813 {
814  // TECs
815  if(beam.isTecInternal(0)){
816  globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
817  globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
818  }
819  // ATs
820  else{
821  // TECminus
822  if(hit < gHitsAtTecMinus){ // gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0
823  globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
824  globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
825  }
826  // TECplus
827  else if(hit > gHitZprime.size() - hitsAtTecPlus - 1){ // gHitZprime[hit] > gBeamSplitterZprime
828  globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhi, globHit[hit].z())));
829  globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z() + 1.0));
830  }
831  // Barrel
832  else{
833  globPtrack.push_back(GlobalPoint(GlobalPoint::Cylindrical(globHit[hit].perp(), trackPhi, globHit[hit].z())));
834  globPref = GlobalPoint(GlobalPoint::Cylindrical(gBeamR, trackPhiRef, globHit[hit].z()));
835  }
836  }
837 }
838 
839 
840 // ----------- create TrajectoryStateOnSurface for each track hit ----------------------------------------------
842  vector<const GeomDetUnit*> &gd, vector<GlobalPoint> &globPtrack,
843  vector<TrajectoryStateOnSurface> &tsosLas, GlobalPoint &globPref)
844 {
846  GlobalVector trajectoryState;
847 
848  // TECplus
849  if(beam.isTecInternal(1)){
850  trajectoryState = GlobalVector(globPref-globPtrack[hit]);
851  }
852  // TECminus
853  else if(beam.isTecInternal(-1)){
854  trajectoryState = GlobalVector(globPtrack[hit]-globPref);
855  }
856  // ATs
857  else{
858  // TECminus
859  if(gHitZprime[hit] < -gBeamSplitterZprime - 2.0*gBeamZ0){
860  trajectoryState = GlobalVector(globPtrack[hit]-globPref);
861  }
862  // TECplus
863  else if(gHitZprime[hit] > gBeamSplitterZprime){
864  trajectoryState = GlobalVector(globPref-globPtrack[hit]);
865  }
866  // Barrel
867  else{
868  trajectoryState = GlobalVector(globPtrack[hit]-globPref);
869  }
870  }
871 // cout << "trajectory: " << trajectoryState << endl;
872  const FreeTrajectoryState ftsLas = FreeTrajectoryState(globPtrack[hit],trajectoryState,0,magneticField);
873  tsosLas.push_back(TrajectoryStateOnSurface(ftsLas,gd[hit]->surface(),
875 }
876 
877 
878 //---------------------- set beam parameters for fittedBeams ---------------------------------
880  unsigned int &hitsAtTecPlus, unsigned int &nFitParams,
881  double &offset, double &slope, vector<GlobalPoint> &globPtrack,
882  double &bsAngleParam, double &chi2)
883 {
884  // set beam parameters for beam output
885  unsigned int paramType(0);
886  if(!fitBeamSplitters_) paramType = 1;
887  if(beam.isAlignmentTube() && hitsAtTecPlus == 0) paramType = 0;
888 // const unsigned int nPedeParams = nFitParams + paramType;
889 
890  // test without BS params
891  const unsigned int nPedeParams(nFitParams);
892 // cout << "number of Pede parameters: " << nPedeParams << endl;
893 
894  std::vector<TkFittedLasBeam::Scalar> params(nPedeParams);
895  params[0] = offset;
896  params[1] = slope;
897  // no BS parameter for AT beams without TECplus hits
898 // if(beam.isTecInternal() || hitsAtTecPlus > 0) params[2] = bsAngleParam;
899 
900  AlgebraicMatrix derivatives(gHitZprime.size(), nPedeParams);
901  // fill derivatives matrix with local track derivatives
902  for(unsigned int hit = 0; hit < gHitZprime.size(); ++hit){
903 
904  // d(delta phi)/d(offset) is identical for every hit
905  derivatives[hit][0] = 1.0;
906 
907  // d(delta phi)/d(slope) and d(delta phi)/d(bsAngleParam) depend on parametrizations
908  // TECplus
909  if(beam.isTecInternal(1)){
910  derivatives[hit][1] = globPtrack[hit].z();
911 // if(gHitZprime[hit] < gBeamSplitterZprime){
912 // derivatives[hit][2] = 0.0;
913 // }
914 // else{
915 // derivatives[hit][2] = - (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
916 // }
917  }
918  // TECminus
919  else if(beam.isTecInternal(-1)){
920  derivatives[hit][1] = globPtrack[hit].z();
921 // if(gHitZprime[hit] > gBeamSplitterZprime){
922 // derivatives[hit][2] = 0.0;
923 // }
924 // else{
925 // derivatives[hit][2] = (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
926 // }
927  }
928  // ATs
929  else{
930  // TECminus
932  derivatives[hit][1] = globPtrack[hit].z();
933 // if(hitsAtTecPlus > 0){
934 // derivatives[hit][2] = 0.0;
935 // }
936  }
937  // TECplus
938  else if(gHitZprime[hit] > gBeamSplitterZprime){
939  derivatives[hit][1] = globPtrack[hit].z();
940 // if(hitsAtTecPlus > 0){
941 // derivatives[hit][2] = - (globPtrack[hit].z() - gBeamSplitterZprime) / gBeamR;
942 // }
943  }
944  // Barrel
945  else{
946  derivatives[hit][1] = globPtrack[hit].z() - gBarrelModuleOffset[hit-gHitsAtTecMinus];
947 // if(hitsAtTecPlus > 0){
948 // derivatives[hit][2] = 0.0;
949 // }
950  }
951  }
952  }
953 
954  unsigned int firstFixedParam(covMatrix.num_col()); // FIXME! --> no, is fine!!!
955 // unsigned int firstFixedParam = nPedeParams - 1;
956 // if(beam.isAlignmentTube() && hitsAtTecPlus == 0) firstFixedParam = nPedeParams;
957 // cout << "first fixed parameter: " << firstFixedParam << endl;
958  // set fit results
959  beam.setParameters(paramType, params, covMatrix, derivatives, firstFixedParam, chi2);
960 
961  return true; // return false in case of problems
962 }
963 
964 
965 //---------------------------------------------------------------------------------------
966 //define this as a plug-in
ESHandle< MagneticField > fieldHandle
static vector< double > gHitZprime
static float gTOBparam
bool getByLabel(std::string const &label, Handle< PROD > &result) const
Definition: Run.h:214
static double gBeamR
static double tecMinusFunction(double *x, double *par)
static double tecPlusFunction(double *x, double *par)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
TkLasBeamFitter(const edm::ParameterSet &config)
static const double slope[3]
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
bool isRing6(void) const
true if this beam hits TEC R6 (last digit of beamId)
Definition: TkLasBeam.h:49
static vector< double > gBarrelModuleRadius
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
virtual LocalPoint localPosition() const
bool fitBeam(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, unsigned int &hitsAtTecPlus, unsigned int &nFitParams, double &offset, double &slope, vector< GlobalPoint > &globPtrack, double &bsAngleParam, double &chi2)
static bool gFitBeamSplitters
bool isTecInternal(int side=0) const
true if this is a TEC internal beam (from 10^2 digit of beamId). side parameter: -1 = ask if TEC-...
Definition: TkLasBeam.cc:4
float float float z
Handle< TkLasBeamCollection > laserBeams
static bool gIsInnerBarrel
int iEvent
Definition: GenABIO.cc:230
virtual void endRunProduce(edm::Run &run, const edm::EventSetup &setup) override
unsigned int nAtParameters_
CLHEP::HepMatrix AlgebraicMatrix
unsigned int getBeamId(void) const
return the full beam identifier
Definition: TkLasBeam.h:25
T z() const
Definition: PV3DBase.h:64
void globalTrackPoint(TkFittedLasBeam &beam, unsigned int &hit, unsigned int &hitsAtTecPlus, double &trackPhi, double &trackPhiRef, std::vector< GlobalPoint > &globHit, std::vector< GlobalPoint > &globPtrack, GlobalPoint &globPref, std::vector< double > &hitPhiError)
static double gBeamSplitterZprime
std::vector< SiStripLaserRecHit2D >::const_iterator end(void) const
access iterator to the collection of hits
Definition: TkLasBeam.h:34
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
virtual void produce(edm::Event &event, const edm::EventSetup &setup) override
std::vector< SiStripLaserRecHit2D >::const_iterator begin(void) const
access iterator to the collection of hits
Definition: TkLasBeam.h:31
How EventSelector::AcceptEvent() decides whether to accept an event for output otherwise it is excluding the probing of A single or multiple positive and the trigger will pass if any such matching triggers are PASS or EXCEPTION[A criterion thatmatches no triggers at all is detected and causes a throw.] A single negative with an expectation of appropriate bit checking in the decision and the trigger will pass if any such matching triggers are FAIL or EXCEPTION A wildcarded negative criterion that matches more than one trigger in the trigger but the state exists so we define the behavior If all triggers are the negative crieriion will lead to accepting the event(this again matches the behavior of"!*"before the partial wildcard feature was incorporated).The per-event"cost"of each negative criterion with multiple relevant triggers is about the same as!*was in the past
static unsigned int gHitsAtTecMinus
void fitter(TkFittedLasBeam &beam, AlgebraicSymMatrix &covMatrix, unsigned int &hitsAtTecPlus, unsigned int &nFitParams, std::vector< double > &hitPhi, std::vector< double > &hitPhiError, std::vector< double > &hitZprimeError, double &zMin, double &zMax, double &bsAngleParam, double &offset, double &offsetError, double &slope, double &slopeError, double &phiAtMinusParam, double &phiAtPlusParam, double &atThetaSplitParam)
const edm::InputTag src_
const SiStripDetId & getDetId(void) const
void buildTrajectory(TkFittedLasBeam &beam, unsigned int &hit, vector< const GeomDetUnit * > &gd, std::vector< GlobalPoint > &globPtrack, vector< TrajectoryStateOnSurface > &tsosLas, GlobalPoint &globPref)
static double gBSparam
const T & get() const
Definition: EventSetup.h:55
T const * product() const
Definition: ESHandle.h:86
static vector< double > gBarrelModuleOffset
ESHandle< TrackerGeometry > geometry
T perp() const
Magnitude of transverse component.
void trackPhi(TkFittedLasBeam &beam, unsigned int &hit, double &trackPhi, double &trackPhiRef, double &offset, double &slope, double &bsAngleParam, double &phiAtMinusParam, double &phiAtPlusParam, double &atThetaSplitParam, std::vector< GlobalPoint > &globHit)
std::vector< TkFittedLasBeam > TkFittedLasBeamCollection
CLHEP::HepSymMatrix AlgebraicSymMatrix
tuple gMinuit
Definition: fitWZ.py:35
static float gTIBparam
static double atFunction(double *x, double *par)
tuple cout
Definition: gather_cfg.py:121
void put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Run.h:107
void getLasBeams(TkFittedLasBeam &beam, vector< TrajectoryStateOnSurface > &tsosLas)
Definition: DDAxes.h:10
std::vector< std::vector< TrajectoryStateOnSurface > > TsosVectorCollection
void getLasHits(TkFittedLasBeam &beam, const SiStripLaserRecHit2D &hit, vector< const GeomDetUnit * > &gd, vector< GlobalPoint > &globHit, unsigned int &hitsAtTecPlus)
int col
Definition: cuy.py:1008
edm::Service< TFileService > fs
std::vector< SiStripLaserRecHit2D >::const_iterator const_iterator
Definition: TkLasBeam.h:16
void setup(std::vector< TH2F > &depth, std::string name, std::string units="")
bool isAlignmentTube(void) const
true if this is an AT beam (from 10^2 digit of beamId)
Definition: TkLasBeam.h:46
Definition: Run.h:41
Global3DVector GlobalVector
Definition: GlobalVector.h:10
list at
Definition: asciidump.py:428
void setParameters(unsigned int parametrisation, const std::vector< Scalar > &params, const AlgebraicSymMatrix &paramCovariance, const AlgebraicMatrix &derivatives, unsigned int firstFixedParam, float chi2)
static double gBeamZ0
Definition: DDAxes.h:10