CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlCaHOCalibProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Feb09 2009
4 // Move the initialisation of SteppingHelixPropagator from ::beginJob() to ::produce()
5 //
6 // Oct3 2008
7 // Difference in tag V00-02-45 with previous code
8 
9 // 1. One new object on data format, which was realised in
10 // CRUZET data analysis.
11 //2. Remove all histogram and cout in the code
12 //3. An upgrade in code, which increases the acceptance of
13 // muon near the edge (this also realised in CRUZET data).
14 // Difference in wrt V00-02-45
15 // 1. initialisation tmpHOCalib.htime = -1000;
16 // 2. By mistake HLT was commented out
17 
18 // Package: AlCaHOCalibProducer
19 // Class: AlCaHOCalibProducer
20 //
53 //
54 // Original Author: Gobinda Majumder
55 // Created: Fri Jul 6 17:17:21 CEST 2007
56 // $Id: AlCaHOCalibProducer.cc,v 1.27 2012/12/26 15:36:24 innocent Exp $
57 //
58 //
59 
60 
61 // system include files
62 #include <memory>
63 
64 // user include files
67 
74 
79 
86 
92 
93 
96 
101 
102 //08/07/07 #include "CondTools/Hcal/interface/HcalDbPool.h"
103 //#include "CondFormats/HcalObjects/interface/HcalPedestals.h"
104 //#include "CondFormats/HcalObjects/interface/HcalPedestalWidths.h"
105 
106 
107 // #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
110 
113 
116 #include "CLHEP/Vector/LorentzVector.h"
117 
120 
125 //#include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixStateInfo.h"
126 
129 
133 
134 
135 #include "TFile.h"
136 #include "TH1F.h"
137 #include "TH2F.h"
138 #include "TProfile.h"
139 #include "TTree.h"
140 /* C++ Headers */
141 #include <string>
142 
143 #include <iostream>
144 #include <fstream>
145 //
146 // class decleration
147 //
148 using namespace std;
149 using namespace edm;
150 using namespace reco;
151 
152 const int netabin= 16;
153 const int nphimx = 72;
154 const int netamx = 32;
155 const int nchnmx = 10;
156 const int ncidmx = 5;
157 
158 //GMA #ifdef DEBUG_OFFLINE_GM
159 const int nsigpk = 7;
160 const int nstrbn = 0;
161 const int ntrgp_gm = 11;
162 
163 
164 const int netahbmx = 60;
165 const int netahb3mx = 32;
166 
167 static const unsigned int nL1trg = 200;
168 
169 static const unsigned int nL1mx=140;
170 static const unsigned int nHLTmx=140;
171 //GMA #endif
172 
174  public:
175  explicit AlCaHOCalibProducer(const edm::ParameterSet&);
177 
181 
182 
183  private:
184  void findHOEtaPhi(int iphsect, int& ietaho, int& iphiho);
185  virtual void beginJob() ;
186  virtual void produce(edm::Event&, const edm::EventSetup&);
187  virtual void endJob() ;
188 
189  // ----------member data ---------------------------
190 
191  float xhor0; //x-position in ring 0
192  float yhor0; //y-position in ring 0
193  float xhor1; //x-position in ring 1
194  float yhor1; //y-position in ring 1
195  int iring; //Ring number -2,-1,0,1,2
196 
197  float localxhor0; //local x-distance from edege in ring 0
198  float localyhor0; //local y-distance from edege in ring 0
199  float localxhor1; //local x-distance from edege in ring 1
200  float localyhor1; //local y-distance from edege in ring 1
201 
202  float pedestal[netamx][nphimx][ncidmx];
203 
205 
206  bool debug;
208 
209  //GMA #ifdef DEBUG_OFFLINE_GM
210 
211  TH1F* libhoped;
212  TH1F* libhoped1;
213 
214  TH1F* allhotime;
215  TH1F* hotime[ntrgp_gm+1];
216  TH1F* hopedtime;
217 
218  TProfile* hopedpr;
219  TH1F* hopedrms;
220  TH1F* hst_hopedrms;
221 
222  TProfile* hopeak[ntrgp_gm+1];
223  TProfile* horatio;
224 
225  TH1F* Nallhotime;
226  TH1F* Nhotime[ntrgp_gm+1];
227  TH1F* Nhopedtime;
228 
229  TH1F* allhb1;
230  TH1F* allhb2;
231  TH1F* allhb3;
232 
233  TH1F* Nallhb1;
234  TH1F* Nallhb2;
235  TH1F* Nallhb3;
236 
237  TProfile* hb1pedpr;
238  TH1F* hb1pedrms;
240 
241  TH1F* ho_occupency[5];
242 
243  bool m_hotime;
244  //GM #endif
245 
246  edm::InputTag muonTags_; // cosmicMuons or standAloneMuons
252 
253  bool m_digiInput; // digi (true) or rechit (false)
254  bool m_hbinfo;
256  int m_endTS;
257  double m_magscale;
258  double m_sigma;
259 
261  //#ifdef DEBUG_OFFLINE_GM
262  // int Nevents;
263  int Noccu;
264  // int Npass;
265  int nRuns;
266  //#endif
267 
268  int irunold;
269  // SteppingHelixPropagator* stepProp;
270  FreeTrajectoryState getFreeTrajectoryState( const reco::Track& tk, const MagneticField* field, int itag, bool dir);
271 
275 
278 
279  unsigned int Ntp; // # of HLT trigger paths (should be the same for all events!)
280  std::map<std::string, bool> fired;
281 
282 };
283 
284 //
285 // constants, enums and typedefs
286 //
287 
288 //
289 // static data member definitions
290 //
291 
292 //
293 // constructors and destructor
294 //
296  : muonTags_(iConfig.getUntrackedParameter<edm::InputTag>("muons"))
297 
298 {
299  //register your products
300 
301  theRootFileName = iConfig.getUntrackedParameter<string>("RootFileName","tmp.root");
302  m_digiInput = iConfig.getUntrackedParameter<bool>("digiInput", true);
303  m_hbinfo = iConfig.getUntrackedParameter<bool>("hbinfo", false);
304  m_startTS = iConfig.getUntrackedParameter<int>("firstTS", 4);
305 
306  m_hotime = iConfig.getUntrackedParameter<bool>("hotime", false);
307 
308  if(m_startTS<0) m_startTS=0;
309  m_endTS = iConfig.getUntrackedParameter<int>("lastTS", 7);
310  if (m_endTS < m_startTS) m_endTS = m_startTS + 3;
311  if (m_endTS >9) m_endTS=9;
312  m_magscale = iConfig.getUntrackedParameter<double>("m_scale", 4.0);
313  m_sigma = iConfig.getUntrackedParameter<double>("sigma", 1.0);
314 
315  hoLabel_ = iConfig.getParameter<edm::InputTag>("hoInput");
316  hbheLabel_ = iConfig.getParameter<edm::InputTag>("hbheInput");
317  hltLabel_ = iConfig.getParameter<edm::InputTag>("hltInput");
318  l1Label_ = iConfig.getParameter<edm::InputTag>("l1Input");
319  towerLabel_ = iConfig.getParameter<edm::InputTag>("towerInput");
320 
321  produces<HOCalibVariableCollection>("HOCalibVariableCollection").setBranchAlias("HOCalibVariableCollection");
322 
323 
324  if (m_hotime) {
326 
327  char title[200];
328  if ( m_digiInput) {
329  libhoped = fs->make<TH1F>("libhoped", "libhoped", ncidmx*netamx*nphimx, -0.5, ncidmx*netamx*nphimx-0.5);
330  libhoped1 = fs->make<TH1F>("libhoped1", "libhoped1", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
331  allhotime = fs->make<TH1F>("allhotime", "allhotime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
332  for (int ij=0; ij<=ntrgp_gm; ij++) {
333  sprintf(title, "hotime_trgp_%i", ij+1);
334  hotime[ij] = fs->make<TH1F>(title, title, nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
335  sprintf(title, "hopeak_trgp_%i", ij+1);
336  hopeak[ij] = fs->make<TProfile>(title, title,netamx*nphimx, -0.5, netamx*nphimx-0.5);
337  }
338 
339  horatio = fs->make<TProfile>("horatio", "horatio",netamx*nphimx, -0.5, netamx*nphimx-0.5);
340  hopedtime = fs->make<TH1F>("hopedtime", "hopedtime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
341 
342  Nallhotime = fs->make<TH1F>("Nallhotime", "Nallhotime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
343  hopedpr = fs->make<TProfile>("hopedpr", "hopedpr", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
344  hopedrms = fs->make<TH1F>("hopedrms", "hopedrms", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
345  hst_hopedrms = fs->make<TH1F>("hst_hopedrms", "hst_hopedrms", 100, 0.0, 0.1);
346  for (int ij=0; ij<=ntrgp_gm; ij++) {
347  sprintf(title, "Nhotime_trgp_%i", ij+1);
348  Nhotime[ij] = fs->make<TH1F>(title, title, nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
349  }
350  Nhopedtime = fs->make<TH1F>("Nhopedtime", "Nhopedtime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
351  allhb1 = fs->make<TH1F>("allhb1", "allhb1", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
352  allhb2 = fs->make<TH1F>("allhb2", "allhb2", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);
353  allhb3 = fs->make<TH1F>("allhb3", "allhb3", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);
354  Nallhb1 = fs->make<TH1F>("Nallhb1", "Nallhb1", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
355  Nallhb2 = fs->make<TH1F>("Nallhb2", "Nallhb2", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);
356  Nallhb3 = fs->make<TH1F>("Nallhb3", "Nallhb3", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);
357  hb1pedpr = fs->make<TProfile>("hb1pedpr", "hb1pedpr", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
358  hb1pedrms = fs->make<TH1F>("hb1pedrms", "hb1pedrms", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
359  hst_hb1pedrms = fs->make<TH1F>("hst_hb1pedrms", "hst_hb1pedrms", 100, 0., 0.1);
360 
361  }
362  for (int i=0; i<5; i++) {
363  sprintf(title, "ho_occupency (>%i #sigma)", i+2);
364  ho_occupency[i] = fs->make<TH1F>(title, title, netamx*nphimx, -0.5, netamx*nphimx-0.5);
365  }
366  }
367 
368 }
369 
371 {
372 
373  // do anything here that needs to be done at desctruction time
374  // (e.g. close files, deallocate resources etc.)
375 
376  if (m_hotime) {
377  // Write the histos to file
378  if ( m_digiInput) {
379  allhotime->Divide(Nallhotime);
380  for (int ij=0; ij<=ntrgp_gm; ij++) {
381  hotime[ij]->Divide(Nhotime[ij]);
382  }
383  hopedtime->Divide(Nhopedtime);
384  libhoped->Scale(1./max(1,nRuns));
385  libhoped1->Scale(1./max(1,nRuns));
386  for (int i=0; i<nchnmx*netamx*nphimx; i++) {
387  float xx = hopedpr->GetBinError(i+1);
388  if (hopedpr->GetBinEntries(i+1) >0) {
389  hopedrms->Fill(i, xx);
390  hst_hopedrms->Fill(xx);
391  }
392  }
393  allhb1->Divide(Nallhb1);
394  allhb2->Divide(Nallhb2);
395  allhb3->Divide(Nallhb3);
396  for (int i=0; i<nchnmx*netahbmx*nphimx; i++) {
397  float xx = hb1pedpr->GetBinError(i+1);
398  if (hb1pedpr->GetBinEntries(i+1) >0) {
399  hb1pedrms->Fill(i, xx);
400  hst_hb1pedrms->Fill(xx);
401  }
402  }
403  }
404  for (int i=0; i<5; i++) {
405  ho_occupency[i]->Scale(1./max(1,Noccu));
406  }
407  }
408 
409 }
410 
411 
412 //
413 // member functions
414 //
415 
416 // ------------ method called to produce the data ------------
417 void
419 {
420 
421  using namespace edm;
422  int irun = iEvent.id().run();
423  if (m_digiInput) {
424  if (irunold !=irun) {
425  iSetup.get<HcalDbRecord>().get(conditions_);
426 
427  for (int i=0; i<netamx; i++) {
428  for (int j=0; j<nphimx; j++) {
429  for (int k=0; k<ncidmx; k++) {
430  pedestal[i][j][k]=0.0;
431  }
432  }
433  }
434  }
435  }
436 
437  // if (m_hotime && m_digiInput) {
438  if (m_digiInput) {
439  if (irunold !=irun) {
440  nRuns++;
441  for (int i =-netabin+1; i<=netabin-1; i++) {
442  if (i==0) continue;
443  int tmpeta1 = (i>0) ? i -1 : -i +14;
444  if (tmpeta1 <0 || tmpeta1 >netamx) continue;
445  for (int j=0; j<nphimx; j++) {
446 
447  HcalDetId id(HcalOuter, i, j+1, 4);
448  calibped = conditions_->getHcalCalibrations(id);
449 
450  for (int k =0; k<ncidmx-1; k++) {
451  pedestal[tmpeta1][j][k] = calibped.pedestal(k); // pedm->getValue(k);
452  pedestal[tmpeta1][j][ncidmx-1] += (1./(ncidmx-1))*pedestal[tmpeta1][j][k];
453  }
454 
455  if (m_hotime) {
456  for (int k =0; k<ncidmx; k++) {
457  libhoped->Fill(nphimx*ncidmx*tmpeta1 + ncidmx*j + k, pedestal[tmpeta1][j][k]);
458  }
459  for (int k =0; k<nchnmx; k++) {
460  libhoped1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*j + k, pedestal[tmpeta1][j][min(k,ncidmx-1)]);
461  }
462  }
463 
464  }
465  }
466  }
467  }
468 
469  // Nevents++;
470  irunold = irun;
471 
472  //GMA if (Nevents%500==1)
473  //GMA cout <<"AlCaHOCalibProducer Processing event # "<<Nevents<<" "<<Npass<<" "<<Noccu<<" "<<irun<<" "<<iEvent.id().event()<<endl;
474 
475  std::auto_ptr<HOCalibVariableCollection> hostore (new HOCalibVariableCollection);
476 
478 
480 
481  if (m_digiInput) {
482  iEvent.getByLabel(hoLabel_,ho);
483  iEvent.getByLabel(hbheLabel_,hbhe);
484  }
485 
486  if (m_hotime && m_digiInput) {
487  if ((*ho).size()>0) {
488  for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
489  HcalDetId id =(*j).id();
490  m_coder = (*conditions_).getHcalCoder(id);
491  m_shape = (*conditions_).getHcalShape(m_coder);
492  int tmpeta= id.ieta();
493  int tmpphi= id.iphi();
494  float tmpdata[nchnmx];
495  int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14;
496  for (int i=0; i<(*j).size() && i<nchnmx; i++) {
497  tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
498  allhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
499  Nallhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.);
500  }
501  }
502  }
503  if ((*hbhe).size()>0) {
504  for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
505  HcalDetId id =(*j).id();
506  m_coder = (*conditions_).getHcalCoder(id);
507  m_shape = (*conditions_).getHcalShape(m_coder);
508  int tmpeta= id.ieta();
509  int tmpphi= id.iphi();
510  int tmpdepth =id.depth();
511  int tmpeta1 = (tmpeta>0) ? tmpeta -15 : -tmpeta + 1;
512  if (tmpdepth==1) tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +29;
513  for (int i=0; i<(*j).size() && i<nchnmx; i++) {
514  float signal = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
515  if (tmpdepth==1) {
516  allhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
517  Nallhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);
518  hb1pedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);}
519  if (tmpdepth==2) {
520  allhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
521  Nallhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
522  if (tmpdepth==3) {
523  allhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
524  Nallhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
525  }
526  }
527  }
528  }
529 
530  double pival = acos(-1.);
531 
533  iEvent.getByLabel(muonTags_, cosmicmuon);
534 
535  if (cosmicmuon->size()>0) {
536 
537  int l1trg = 0;
538  int hlttr = 0;
539 
540  int ntrgpas_gm[ntrgp_gm]={0,0,0,0,0,0,0,0,0,0};
541 
542  /*
543  //L1 trigger
544  Handle<L1GlobalTriggerReadoutRecord> L1GTRR;
545  iEvent.getByLabel(l1Label_,L1GTRR); //gtDigis
546 
547  if ( L1GTRR.isValid()) {
548  const unsigned int n(L1GTRR->decisionWord().size());
549  const bool accept(L1GTRR->decision());
550  if (accept) {
551  for (unsigned int i=0; i!=n && i<32; ++i) {
552  // for (unsigned int i=0; i!=n ; ++i) {
553  int il1trg = (L1GTRR->decisionWord()[i]) ? 1 : 0;
554  if (il1trg>0 && i<32) l1trg +=int(std::pow(2., double(i%32))*il1trg);
555  }
556  }
557  }// else { return;}
558 
559  //HLT
560 
561  Handle<edm::TriggerResults> trigRes;
562  iEvent.getByLabel(hltLabel_, trigRes);
563 
564 
565  unsigned int size = trigRes->size();
566  edm::TriggerNames triggerNames(*trigRes);
567 
568  // loop over all paths, get trigger decision
569  for(unsigned i = 0; i != size && i<32; ++i) {
570  std::string name = triggerNames.triggerName(i);
571  fired[name] = trigRes->accept(i);
572  int ihlt = trigRes->accept(i);
573  if (m_hotime){
574  if (ihlt >0 && i < (int)ntrgp_gm) { ntrgpas_gm[i] = 1;}
575  }
576  if (i<32 && ihlt>0) hlttr += int(std::pow(2., double(i%32))*ihlt);
577  }
578 
579  */
580 
581  int Noccu_old = Noccu;
582 
583  for(reco::TrackCollection::const_iterator ncosm = cosmicmuon->begin();
584  ncosm != cosmicmuon->end(); ++ncosm) {
585 
586  if ((*ncosm).ndof() < 15) continue;
587  if ((*ncosm).normalizedChi2() >30.0) continue;
588 
589  HOCalibVariables tmpHOCalib;
590 
591  tmpHOCalib.trig1 = l1trg;
592  tmpHOCalib.trig2 = hlttr;
593 
594  int charge = ncosm->charge();
595 
596  double innerr = (*ncosm).innerPosition().Perp2();
597  double outerr = (*ncosm).outerPosition().Perp2();
598  int iiner = (innerr <outerr) ? 1 : 0;
599 
600  //---------------------------------------------------
601  // in_to_out Dir in_to_out Dir
602  // StandAlone ^ ^ Cosmic ^ |
603  // | | | v
604  //---------------------------------------------------Y=0
605  // StandAlone | | Cosmic ^ |
606  // v v | v
607  //----------------------------------------------------
608 
609  double posx, posy, posz;
610  double momx, momy, momz;
611 
612  if (iiner==1) {
613  posx = (*ncosm).innerPosition().X();
614  posy = (*ncosm).innerPosition().Y();
615  posz = (*ncosm).innerPosition().Z();
616 
617  momx = (*ncosm).innerMomentum().X();
618  momy = (*ncosm).innerMomentum().Y();
619  momz = (*ncosm).innerMomentum().Z();
620 
621  } else {
622  posx = (*ncosm).outerPosition().X();
623  posy = (*ncosm).outerPosition().Y();
624  posz = (*ncosm).outerPosition().Z();
625 
626  momx = (*ncosm).outerMomentum().X();
627  momy = (*ncosm).outerMomentum().Y();
628  momz = (*ncosm).outerMomentum().Z();
629  }
630 
631 
632  PositionType trkpos(posx, posy, posz);
633 
634  CLHEP::Hep3Vector tmpmuon3v(posx, posy, posz);
635  CLHEP::Hep3Vector tmpmuondir(momx, momy, momz);
636 
637  bool samedir = (tmpmuon3v.dot(tmpmuondir) >0) ? true : false;
638  for (int i=0; i<3; i++) {tmpHOCalib.caloen[i] = 0.0;}
639  int inearbymuon = 0;
640  for(reco::TrackCollection::const_iterator ncosmcor = cosmicmuon->begin();
641  ncosmcor != cosmicmuon->end(); ++ncosmcor) {
642  if (ncosmcor==ncosm) continue;
643 
644  CLHEP::Hep3Vector tmpmuon3vcor;
645  CLHEP::Hep3Vector tmpmom3v;
646  if (iiner==1) {
647  tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).innerPosition().X(),(*ncosmcor).innerPosition().Y(),(*ncosmcor).innerPosition().Z());
648  tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).innerMomentum().X(),(*ncosmcor).innerMomentum().Y(),(*ncosmcor).innerMomentum().Z());
649  } else {
650  tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).outerPosition().X(),(*ncosmcor).outerPosition().Y(),(*ncosmcor).outerPosition().Z());
651  tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).outerMomentum().X(),(*ncosmcor).outerMomentum().Y(),(*ncosmcor).outerMomentum().Z());
652 
653  }
654  if (tmpmom3v.mag()<0.2 || (*ncosmcor).ndof()<5) continue;
655 
656  double angle = tmpmuon3v.angle(tmpmuon3vcor);
657  if (angle < 7.5*pival/180.) {inearbymuon=1;} // break;}
658 
659  if (muonTags_.label() =="cosmicMuons") {
660  if (angle <7.5*pival/180.) { tmpHOCalib.caloen[0] +=1.;}
661  if (angle <15.0*pival/180.) { tmpHOCalib.caloen[1] +=1.;}
662  if (angle <35.0*pival/180.) { tmpHOCalib.caloen[2] +=1.;}
663  }
664  }
665 
666  localxhor0 = localyhor0 = 20000; //GM for 22OCT07 data
667 
668  if (muonTags_.label() =="standAloneMuons") {
669 
670  Handle<CaloTowerCollection> calotower;
671  iEvent.getByLabel(towerLabel_, calotower);
672 
673  for (CaloTowerCollection::const_iterator calt = calotower->begin();
674  calt !=calotower->end(); calt++) {
675  //CMSSW_2_1_x const math::XYZVector towermom = (*calt).momentum();
676  double ith = (*calt).momentum().theta();
677  double iph = (*calt).momentum().phi();
678 
679  CLHEP::Hep3Vector calo3v(sin(ith)*cos(iph), sin(ith)*sin(iph), cos(ith));
680 
681  double angle = tmpmuon3v.angle(calo3v);
682 
683  if (angle < 7.5*pival/180.) {tmpHOCalib.caloen[0] += calt->emEnergy()+calt->hadEnergy();}
684  if (angle < 15*pival/180.) {tmpHOCalib.caloen[1] += calt->emEnergy()+calt->hadEnergy();}
685  if (angle < 35*pival/180.) {tmpHOCalib.caloen[2] += calt->emEnergy()+calt->hadEnergy();}
686  }
687 
688 
689  }
690  if (tmpHOCalib.caloen[0] >10.0) continue;
691 
692  GlobalPoint glbpt(posx, posy, posz);
693 
694  double mom = sqrt(momx*momx + momy*momy +momz*momz);
695 
696  momx /= mom;
697  momy /= mom;
698  momz /= mom;
699 
700  DirectionType trkdir(momx, momy, momz);
701 
702  tmpHOCalib.trkdr = (*ncosm).d0();
703  tmpHOCalib.trkdz = (*ncosm).dz();
704 
705  tmpHOCalib.nmuon = cosmicmuon->size();
706  tmpHOCalib.trkvx = glbpt.x();
707  tmpHOCalib.trkvy = glbpt.y();
708  tmpHOCalib.trkvz = glbpt.z();
709  tmpHOCalib.trkmm = mom*charge;
710  tmpHOCalib.trkth = trkdir.theta();
711  tmpHOCalib.trkph = trkdir.phi();
712 
713  tmpHOCalib.ndof = (inearbymuon ==0) ? (int)(*ncosm).ndof() : -(int)(*ncosm).ndof();
714  tmpHOCalib.chisq = (*ncosm).normalizedChi2(); // max(1.,tmpHOCalib.ndof);
715  tmpHOCalib.therr = 0.;
716  tmpHOCalib.pherr = 0.;
717 
718  if (iiner==1) {
719  reco::TrackBase::CovarianceMatrix innercov = (*ncosm).innerStateCovariance();
720  tmpHOCalib.therr = innercov(1,1); //thetaError();
721  tmpHOCalib.pherr = innercov(2,2); //phi0Error();
722  } else {
723  reco::TrackBase::CovarianceMatrix outercov = (*ncosm).outerStateCovariance();
724  tmpHOCalib.therr = outercov(1,1); //thetaError();
725  tmpHOCalib.pherr = outercov(2,2); //phi0Error();
726  }
727 
728  ESHandle<MagneticField> theMagField;
729  iSetup.get<IdealMagneticFieldRecord>().get(theMagField );
730 
731  SteppingHelixPropagator myHelix(&*theMagField,anyDirection);
732  myHelix.setMaterialMode(false);
733  myHelix.applyRadX0Correction(true);
734 
735  double phiho = trkpos.phi();
736  if (phiho<0) phiho +=2*pival;
737 
738  int iphisect_dt=int(6*(phiho+pival/18.)/pival); //for u 18/12/06
739  if (iphisect_dt>=12) iphisect_dt=0;
740 
741  int iphisect = -1;
742 
743  bool ipath = false;
744  for (int kl = 0; kl<=2; kl++) {
745 
746  int iphisecttmp = (kl<2) ? iphisect_dt + kl : iphisect_dt - 1;
747  if (iphisecttmp <0) iphisecttmp = 11;
748  if (iphisecttmp >=12) iphisecttmp = 0;
749 
750  double phipos = iphisecttmp*pival/6.;
751  double phirot = phipos;
752 
753  GlobalVector xLocal(-sin(phirot), cos(phirot), 0.);
754 
755  GlobalVector yLocal(0., 0., 1.);
756  GlobalVector zLocal = xLocal.cross(yLocal).unit();
757  // GlobalVector zLocal(cos(phirot), sin(phirot), 0.0);
758 
759 
760  FreeTrajectoryState freetrajectorystate_ = getFreeTrajectoryState(*ncosm,&(*theMagField), iiner, samedir);
761 
762  Surface::RotationType rot(xLocal, yLocal, zLocal);
763 
764  for (int ik=1; ik>=0; ik--) { //propagate track in two HO layers
765 
766  double radial = 407.0;
767  if (ik==0) radial = 382.0;
768 
769  Surface::PositionType pos(radial*cos(phipos), radial*sin(phipos), 0.);
770  PlaneBuilder::ReturnType aPlane = PlaneBuilder().plane(pos,rot);
771 
772  auto aPlane2 = new Plane(pos,rot);
773 
774  SteppingHelixStateInfo steppingHelixstateinfo_ = myHelix.propagate(SteppingHelixStateInfo(freetrajectorystate_), (*aPlane2));
775 
776  if (steppingHelixstateinfo_.isValid()) {
777 
778  GlobalVector hotrkpos2(steppingHelixstateinfo_.position().x(), steppingHelixstateinfo_.position().y(), steppingHelixstateinfo_.position().z());
779  CLHEP::Hep3Vector hotrkdir2(steppingHelixstateinfo_.momentum().x(), steppingHelixstateinfo_.momentum().y(),steppingHelixstateinfo_.momentum().z());
780 
781  LocalVector lclvt0 = (*aPlane).toLocal(hotrkpos2);
782 
783  double xx = lclvt0.x();
784  double yy = lclvt0.y();
785 
786  if (ik ==1) {
787  if ((std::abs(yy) < 130 && xx >-64.7 && xx <138.2)
788  ||(std::abs(yy) > 130 && std::abs(yy) <700 && xx >-76.3 && xx <140.5)) {
789  ipath = true; //Only look for tracks which as hits in layer 1
790  iphisect = iphisecttmp;
791  }
792  }
793 
794  if (iphisect != iphisecttmp) continue; //Look for ring-0 only when ring1 is accepted for that sector
795 
796  switch (ik)
797  {
798  case 0 :
799  xhor0 = xx; //lclvt0.x();
800  yhor0 = yy; //lclvt0.y();
801  break;
802  case 1 :
803  xhor1 = xx; //lclvt0.x();
804  yhor1 = yy; //lclvt0.y();
805 
806  tmpHOCalib.hoang = CLHEP::Hep3Vector(zLocal.x(),zLocal.y(),zLocal.z()).dot(hotrkdir2.unit());
807  break;
808  default : break;
809  }
810  } else {
811  break;
812  }
813  }
814  if (ipath) break;
815  }
816  if (ipath) { //If muon crossed HO laeyrs
817 
818  int ietaho = 50;
819  int iphiho = -1;
820 
821  for (int i=0; i<9; i++) {tmpHOCalib.hosig[i]=-100.0;}
822  for (int i=0; i<18; i++) {tmpHOCalib.hocorsig[i]=-100.0;}
823  for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
824  tmpHOCalib.hocro = -100;
825  tmpHOCalib.htime = -1000;
826 
827  int isect = 0;
828 
829  findHOEtaPhi(iphisect, ietaho, iphiho);
830 
831  if (ietaho !=0 && iphiho !=0 && std::abs(iring)<=2) { //Muon passed through a tower
832  isect = 100*std::abs(ietaho+30)+std::abs(iphiho);
833  if (std::abs(ietaho) >=netabin || iphiho<0) isect *=-1; //Not extrapolated to any tower
834  if (std::abs(ietaho) >=netabin) isect -=1000000; //not matched with eta
835  if (iphiho<0) isect -=2000000; //not matched with phi
836  tmpHOCalib.isect = isect;
837 
838  tmpHOCalib.hodx = localxhor1;
839  tmpHOCalib.hody = localyhor1;
840 
841  if (iring==0) {
842  tmpHOCalib.hocorsig[8] = localxhor0;
843  tmpHOCalib.hocorsig[9] = localyhor0;
844  }
845 
846  int etamn=-4;
847  int etamx=4;
848  if (iring==1) {etamn=5; etamx = 10;}
849  if (iring==2) {etamn=11; etamx = 16;}
850  if (iring==-1){etamn=-10; etamx = -5;}
851  if (iring==-2){etamn=-16; etamx = -11;}
852 
853  int phimn = 1;
854  int phimx = 2;
855  if (iring ==0) {
856  phimx =2*int((iphiho+1)/2.);
857  phimn = phimx - 1;
858  } else {
859  phimn = 3*int((iphiho+1)/3.) - 1;
860  phimx = phimn + 2;
861  }
862 
863  if (phimn <1) phimn += nphimx;
864  if (phimx >72) phimx -= nphimx;
865 
866  int sigstr = m_startTS; // 5;
867  int sigend = m_endTS; // 8;
868 
869  // if (iphiho <=nphimx/2) { //GMA310508
870  // sigstr -=1; //GMA300608
871  // sigend -=1;
872  // }
873 
874  if (m_hbinfo) {
875  for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
876 
877  if (m_digiInput) {
878  if ((*hbhe).size() >0) {
879  for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
880  // const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
881  // HcalDetId id =digi.id();
882  HcalDetId id =(*j).id();
883  m_coder = (*conditions_).getHcalCoder(id);
884  m_shape = (*conditions_).getHcalShape(m_coder);
885  int tmpeta= id.ieta();
886  int tmpphi= id.iphi();
887  calibped = conditions_->getHcalCalibrations(id);
888 
889  int deta = tmpeta-ietaho;
890  if (tmpeta==-1 && ietaho== 1) deta = -1;
891  if (tmpeta== 1 && ietaho==-1) deta = 1;
892  int dphi = tmpphi-iphiho;
893  if (phimn >phimx) {
894  if (dphi==71) dphi=-1;
895  if (dphi==-71) dphi=1;
896  }
897 
898  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
899 
900  if (ipass2 ==0 ) continue;
901 
902  float tmpdata[nchnmx];
903  for (int i=0; i<(*j).size() && i<nchnmx; i++) {
904  tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
905  }
906 
907  float signal = 0;
908  for (int i=1; i<(*j).size() && i<=8; i++) {
909  signal += tmpdata[i] - calibped.pedestal((*j).sample(i).capid());;
910  }
911 
912  if (ipass2 == 1) {
913  if (3*(deta+1)+dphi+1<9) tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
914  }
915  }
916  }
917 
918  } else {
919 
920  edm::Handle<HBHERecHitCollection> hbheht;// iEvent.getByType(hbheht);
921  iEvent.getByLabel(hbheLabel_,hbheht);
922 
923 
924  if ((*hbheht).size()>0) {
925  if(!(*hbheht).size()) throw (int)(*hbheht).size();
926 
927  for (HBHERecHitCollection::const_iterator j=(*hbheht).begin(); j!=(*hbheht).end(); j++){
928  // const HBHERecHit hbhehtrec = (const HBHERecHit)(*j);
929  // HcalDetId id =hbhehtrec.id();
930  HcalDetId id =(*j).id();
931  int tmpeta= id.ieta();
932  int tmpphi= id.iphi();
933 
934  int deta = tmpeta-ietaho;
935  if (tmpeta==-1 && ietaho== 1) deta = -1;
936  if (tmpeta== 1 && ietaho==-1) deta = 1;
937  int dphi = tmpphi-iphiho;
938  if (phimn >phimx) {
939  if (dphi==71) dphi=-1;
940  if (dphi==-71) dphi=1;
941  }
942 
943  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
944  if ( ipass2 ==0 ) continue;
945 
946  float signal = (*j).energy();
947 
948  if (3*(deta+1)+dphi+1<9) tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
949  }
950  }
951 
952  } //else m_digilevel
953 
954  } //m_hbinfo #endif
955 
956  if (m_digiInput) {
957  if ((*ho).size()>0) {
958  int isFilled[netamx*nphimx];
959  for (int j=0; j<netamx*nphimx; j++) {isFilled[j]=0;}
960 
961  double sumEt = 0;
962  double sumE = 0;
963 
964  for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
965  // const HODataFrame digi = (const HODataFrame)(*j);
966  // HcalDetId id =digi.id();
967 
968  HcalDetId id =(*j).id();
969  m_coder = (*conditions_).getHcalCoder(id);
970  m_shape = (*conditions_).getHcalShape(m_coder);
971  int tmpeta= id.ieta();
972  int tmpphi= id.iphi();
973 
974  int ipass1 =0;
975  if (tmpeta >=etamn && tmpeta <=etamx) {
976  if (phimn < phimx) {
977  ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
978  } else {
979  ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
980  }
981  }
982 
983  int deta = tmpeta-ietaho;
984  if (tmpeta==-1 && ietaho== 1) deta = -1;
985  if (tmpeta== 1 && ietaho==-1) deta = 1;
986 
987  int dphi = tmpphi -iphiho;
988  if (phimn>phimx) {
989  if (dphi==71) dphi=-1;
990  if (dphi==-71) dphi=1;
991  }
992 
993  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0;
994 
995  int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14;
996 
997  float tmpdata[nchnmx]={0,0,0,0,0,0,0,0,0,0};
998  float sigvall[nsigpk]={0,0,0,0,0,0,0};
999 
1000  for (int i=0; i<(*j).size() && i<nchnmx; i++) {
1001  tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
1002  if (deta==0 && dphi==0) {
1003  double tmpE = tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
1004  if (tmpE >0) {
1005  sumEt += i*tmpE;
1006  sumE += tmpE;
1007  }
1008  if (m_hotime) {
1009  //calculate signals in 4 time slices, 0-3,.. 6-9
1010  if (i>=7-nsigpk) {
1011  for (int ncap=0; ncap<nsigpk; ncap++) {
1012  if (i-ncap >= nstrbn && i-ncap <= nstrbn+3) {
1013  sigvall[ncap] +=tmpdata[i];
1014  }
1015  }
1016  }
1017  if (i==(*j).size()-1) {
1018  float mxled=-1;
1019  int imxled = 0;
1020  for (int ij=0; ij<nsigpk; ij++) {
1021  if (sigvall[ij] > mxled) {mxled = sigvall[ij]; imxled=ij;}
1022  }
1023  double pedx = 0.0;
1024  for (int ij=0; ij<4; ij++) {
1025  pedx +=pedestal[tmpeta1][tmpphi-1][ij];
1026  }
1027  if (mxled-pedx >2 && mxled-pedx <20 ) {
1028  hopeak[ntrgp_gm]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
1029  for (int jk=0; jk<ntrgp_gm; jk++) {
1030  if (ntrgpas_gm[jk]>0) {
1031  hopeak[jk]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
1032  }
1033  }
1034  if (tmpdata[5]+tmpdata[6] >1) {
1035  horatio->Fill(nphimx*tmpeta1 + tmpphi-1, (tmpdata[5]-tmpdata[6])/(tmpdata[5]+tmpdata[6]));
1036  }
1037  for (int ij=0; ij<(*j).size() && ij<nchnmx; ij++) {
1038  hotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
1039  Nhotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
1040  for (int jk=0; jk<ntrgp_gm; jk++) {
1041  if (ntrgpas_gm[jk]>0) {
1042  hotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
1043  Nhotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
1044  }
1045  }
1046  }
1047  }
1048  }
1049  }
1050  }
1051  }
1052 
1053  if (std::abs(tmpeta) <=15 && deta==0 && dphi ==0) {
1054  float signal = 0;
1055  size_t icnt = 0;
1056  for (int i =0; i<nchnmx && i< (*j).size(); i++) {
1057  if (i >=sigstr && i<=sigend) continue;
1058  signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
1059  if (++icnt >=4) break;
1060  }
1061  tmpHOCalib.hocro = signal;
1062  }
1063 
1064  if (m_hotime) {
1065  if (ipass1 ==0 && ipass2 ==0 && cosmicmuon->size()<=2) {
1066  if (std::abs(ietaho) <=netabin && iphiho >0) {
1067  if ((iphiho >=1 && iphiho<=nphimx/2 && tmpphi >=1 && tmpphi <=nphimx/2) ||
1068  (iphiho >nphimx/2 && iphiho<=nphimx && tmpphi >nphimx/2 && tmpphi <=nphimx)) {
1069  if (isFilled[nphimx*tmpeta1+tmpphi-1]==0) {
1070  isFilled[nphimx*tmpeta1+tmpphi-1]=1;
1071  for (int i=0; i<(*j).size() && i<nchnmx; i++) {
1072  hopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
1073  Nhopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.);
1074  hopedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
1075  }
1076  } //isFilled
1077  }
1078  }
1079  }
1080  }
1081 
1082  if (ipass1 ==0 && ipass2 ==0 ) continue;
1083 
1084  float signal = 0;
1085  for (int i=sigstr; i<(*j).size() && i<=sigend; i++) {
1086  signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
1087  }
1088  if (signal <-100 || signal >100000) signal = -100;
1089  if (m_hotime) {
1090  if (signal >-100 && Noccu == Noccu_old) {
1091  for (int i=0; i<5; i++) {
1092  if (signal >(i+2)*m_sigma) {
1093  ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
1094  }
1095  }
1096  }
1097  }
1098 
1099  if (ipass1 ==0 && ipass2 ==0 ) continue;
1100 
1101  if (ipass1 ==1) {
1102  int tmpdph = tmpphi-phimn;
1103  if (tmpdph<0) tmpdph = 2; //only case of iphi==1, where phimn=71
1104 
1105  int ilog = 2*(tmpeta-etamn)+tmpdph;
1106  if (iring !=0) {
1107  if (iring >0) {
1108  ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
1109  } else {
1110  ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
1111  }
1112  }
1113  if (ilog>-1 && ilog<18) {
1114  tmpHOCalib.hocorsig[ilog] = signal;
1115  }
1116  }
1117 
1118  if (ipass2 ==1) {
1119  if (3*(deta+1)+dphi+1<9) tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
1120  }
1121 
1122  /*
1123  // Possibility to store pedestal by shifting phi tower by 6
1124  // But, due to missing tower at +-5, we do not have always proper
1125  // statistics and also in pedestal subtracted data, we do not have
1126  // signal in that tower
1127  //
1128  if (deta==0 && dphi ==0) {
1129  int crphi = tmpphi + 6;
1130  if (crphi >72) crphi -=72;
1131 
1132  for (HODigiCollection::const_iterator jcr=(*ho).begin(); jcr!=(*ho).end(); jcr++){
1133  // const HODataFrame (*jcr) = (const HODataFrame)(*jcr);
1134  // HcalDetId idcr =(*jcr).id();
1135  HcalDetId id =(*jcr).id();
1136  m_coder = (*conditions_).getHcalCoder(idcr);
1137  m_shape = (*conditions_).getHcalShape(m_coder);
1138  int etacr= idcr.ieta();
1139  int phicr= idcr.iphi();
1140 
1141  if (tmpeta==etacr && crphi ==phicr) {
1142 
1143  float tmpdatacr[nchnmx];
1144  for (int i=0; i<(*jcr).size() && i<nchnmx; i++) {
1145  tmpdatacr[i] = m_coder->charge(*m_shape,(*jcr).sample(i).adc(),(*jcr).sample(i).capid());
1146  }
1147  }
1148  }
1149  }
1150  */
1151 
1152  }
1153  tmpHOCalib.htime = sumEt/max(sumE,1.e-6);
1154  }
1155  } else {
1157  iEvent.getByLabel(hoLabel_,hoht);
1158 
1159 
1160  if ((*hoht).size()>0) {
1161  for (HORecHitCollection::const_iterator j=(*hoht).begin(); j!=(*hoht).end(); j++){
1162  // const HORecHit hohtrec = (const HORecHit)(*j);
1163  // HcalDetId id =hohtrec.id();
1164  HcalDetId id =(*j).id();
1165  int tmpeta= id.ieta();
1166  int tmpphi= id.iphi();
1167 
1168  int ipass1 =0;
1169  if (tmpeta >=etamn && tmpeta <=etamx) {
1170  if (phimn < phimx) {
1171  ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
1172  } else {
1173  ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
1174  }
1175  }
1176 
1177  int deta = tmpeta-ietaho;
1178  if (tmpeta==-1 && ietaho== 1) deta = -1;
1179  if (tmpeta== 1 && ietaho==-1) deta = 1;
1180 
1181  int dphi = tmpphi -iphiho;
1182  if (phimn>phimx) {
1183  if (dphi==71) dphi=-1;
1184  if (dphi==-71) dphi=1;
1185  }
1186 
1187  float signal = (*j).energy();
1188  if (m_hotime) {
1189  int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14;
1190  if (signal >-100 && Noccu == Noccu_old) {
1191  for (int i=0; i<5; i++) {
1192  if (signal >(i+2)*m_sigma) {
1193  ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
1194  }
1195  }
1196  }
1197  }
1198 
1199  int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0;
1200 
1201  if (ipass1 ==0 && ipass2 ==0 ) continue;
1202 
1203  if (ipass1 ==1) {
1204  int tmpdph = tmpphi-phimn;
1205  if (tmpdph<0) tmpdph = 2; //only case of iphi==1, where phimn=71
1206 
1207  int ilog = 2*(tmpeta-etamn)+tmpdph;
1208  if (iring !=0) {
1209  if (iring >0) {
1210  ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
1211  } else {
1212  ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
1213  }
1214  }
1215  if (ilog>-1 && ilog<18) {
1216  tmpHOCalib.hocorsig[ilog] = signal;
1217  }
1218  }
1219 
1220  if (ipass2 ==1) {
1221 
1222  if (3*(deta+1)+dphi+1<9) {
1223  tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
1224  }
1225  }
1226 
1227  if (deta==0 && dphi ==0) {
1228  tmpHOCalib.htime = (*j).time();
1229  int crphi = tmpphi + 6;
1230  if (crphi >72) crphi -=72;
1231 
1232  for (HORecHitCollection::const_iterator jcr=(*hoht).begin(); jcr!=(*hoht).end(); jcr++){
1233  const HORecHit reccr = (const HORecHit)(*jcr);
1234  HcalDetId idcr =reccr.id();
1235  int etacr= idcr.ieta();
1236  int phicr= idcr.iphi();
1237  if (tmpeta==etacr && crphi ==phicr) {
1238 
1239  tmpHOCalib.hocro = reccr.energy();
1240 
1241  }
1242  }
1243  }
1244  }
1245  }
1246  }
1247 
1248  //GMA Npass++;
1249  if (Noccu == Noccu_old) Noccu++;
1250  hostore->push_back(tmpHOCalib);
1251 
1252  }
1253  }
1254 
1255  }
1256  }
1257 
1258  iEvent.put(hostore, "HOCalibVariableCollection");
1259 
1260 }
1261 
1262 // ------------ method called once each job just before starting event loop ------------
1263 void
1265 {
1266  //GMA Nevents = 0;
1267  //GMA Npass = 0;
1268  //GMA Noccu = 0;
1269 
1270  irunold = -1;
1271  nRuns = 0;
1272  // edm::ESHandle<MagneticField> bField;
1273  // iSetup.get<IdealMagneticFieldRecord>().get(bField);
1274  // stepProp = new SteppingHelixPropagator(&*bField,anyDirection);
1275  // stepProp->setMaterialMode(false);
1276  // stepProp->applyRadX0Correction(true);
1277 
1278  for (int i=0; i<netamx; i++) {
1279  for (int j=0; j<nphimx; j++) {
1280  for (int k=0; k<ncidmx; k++) {
1281  pedestal[i][j][k]=0.0;
1282  }
1283  }
1284  }
1285 
1286 
1287 }
1288 
1289 // ------------ method called once each job just after ending the event loop ------------
1290 void
1292 
1293 
1294 }
1295 
1296 void AlCaHOCalibProducer::findHOEtaPhi(int iphisect, int& ietaho, int& iphiho) {
1297 
1298  //18/12/06 : use only position, not angle phi
1299 
1300 double etalow[netabin]={ 0.025, 35.195, 70.625, 106.595, 141.565, 180.765, 220.235, 261.385, 304.525, 349.975, 410.025, 452.085, 506.645, 565.025, 627.725, 660.25};
1301 double etahgh[netabin]={ 35.145, 70.575, 106.545, 125.505, 180.715, 220.185, 261.335, 304.475, 349.925, 392.575, 452.035, 506.595, 564.975, 627.675, 661.075, 700.25};
1302 
1303  double philow[6]={-76.27, -35.11, 0.35, 35.81, 71.77, 108.93}; //Ring+/-1 & 2
1304  double phihgh[6]={-35.81, -0.35, 35.11, 71.07, 108.23, 140.49};
1305 
1306  double philow00[6]={-60.27, -32.91, 0.35, 33.61, 67.37, 102.23}; //Ring0 L0
1307  double phihgh00[6]={-33.61, -0.35, 32.91, 66.67, 101.53, 129.49};
1308 
1309  double philow01[6]={-64.67, -34.91, 0.35, 35.61, 71.37, 108.33}; //Ring0 L1
1310  double phihgh01[6]={-35.61, -0.35, 34.91, 70.67, 107.63, 138.19};
1311 
1312 
1313  iring = -10;
1314 
1315  double tmpdy = std::abs(yhor1);
1316  for (int i=0; i<netabin; i++) {
1317  if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
1318  ietaho = i+1;
1319  float tmp1 = fabs(tmpdy-etalow[i]);
1320  float tmp2 = fabs(tmpdy-etahgh[i]);
1321 
1322  localyhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
1323 
1324  if (i<4) iring =0;
1325  if (i>=4 && i<10) iring=1;
1326  if (i>=10 && i<netabin) iring=2;
1327  break;
1328  }
1329  }
1330 
1331  int tmpphi = 0;
1332  int tmpphi0 = 0;
1333 
1334  if (ietaho >4) { //Ring 1 and 2
1335  for (int i=0; i<6; i++) {
1336  if (xhor1 >philow[i] && xhor1 <phihgh[i]) {
1337  tmpphi=i+1;
1338  float tmp1 = fabs(xhor1-philow[i]);
1339  float tmp2 = fabs(xhor1-phihgh[i]);
1340  localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
1341  break;
1342  }
1343  }
1344  } else { //Ring 0
1345  for (int i=0; i<6; i++) {
1346  if (xhor1 >philow01[i] && xhor1 <phihgh01[i]) {
1347  tmpphi=i+1;
1348  float tmp1 = fabs(xhor1-philow01[i]);
1349  float tmp2 = fabs(xhor1-phihgh01[i]);
1350  localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
1351  break;
1352  }
1353  }
1354 
1355  for (int i=0; i<6; i++) {
1356  if (xhor0 >philow00[i] && xhor0 <phihgh00[i]) {
1357  tmpphi0=i+1;
1358  float tmp1 = fabs(xhor0-philow00[i]);
1359  float tmp2 = fabs(xhor0-phihgh00[i]);
1360  localxhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
1361  if (tmpphi !=tmpphi0) localxhor0 +=10000.;
1362  break;
1363  }
1364  }
1365 
1366  double tmpdy = std::abs(yhor0);
1367  for (int i=0; i<4; i++) {
1368  if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
1369  float tmp1 = fabs(tmpdy-etalow[i]);
1370  float tmp2 = fabs(tmpdy-etahgh[i]);
1371  localyhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
1372  if (i+1 != ietaho) localyhor0 +=10000.;
1373  break;
1374  }
1375  }
1376  }
1377 
1378  if (tmpphi!=0) {
1379  iphiho = 6*iphisect -2 + tmpphi;
1380  if (iphiho <=0) iphiho +=nphimx;
1381  if (iphiho >nphimx) iphiho -=nphimx;
1382  }
1383 
1384  // isect2 = 15*iring+iphisect+1;
1385 
1386  if (yhor1 <0) {
1387  if (std::abs(ietaho) >netabin) { //Initialised with 50
1388  ietaho +=1;
1389  } else {
1390  ietaho *=-1;
1391  }
1392  // isect2 *=-1;
1393  iring *=-1;
1394  }
1395 }
1396 
1398 {
1399 
1400  if (iiner ==0) {
1401  GlobalPoint gpos( tk.outerX(), tk.outerY(), tk.outerZ());
1402  GlobalVector gmom( tk.outerPx(), tk.outerPy(), tk.outerPz());
1403  if (dir) gmom *=-1.;
1404  GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
1405  CurvilinearTrajectoryError err( tk.extra()->outerStateCovariance());
1406  return FreeTrajectoryState( par, err);
1407  } else {
1408  GlobalPoint gpos( tk.innerPosition().X(), tk.innerPosition().Y(), tk.innerPosition().Z());
1409  GlobalVector gmom( tk.innerMomentum().X(), tk.innerMomentum().Y(), tk.innerMomentum().Z());
1410  if (dir) gmom *=-1.;
1411  GlobalTrajectoryParameters par( gpos, -gmom, tk.charge(), field);
1412  CurvilinearTrajectoryError err( tk.extra()->innerStateCovariance());
1413  return FreeTrajectoryState( par, err);
1414  }
1415 
1416 }
1417 
1419 
1420 //define this as a plug-in
1422 
RunNumber_t run() const
Definition: EventID.h:42
T getParameter(std::string const &) const
Basic3DVector< float > DirectionType
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
TProfile * hopeak[ntrgp_gm+1]
FreeTrajectoryState getFreeTrajectoryState(const reco::Track &tk, const MagneticField *field, int itag, bool dir)
const int netabin
double outerPy() const
y coordinate of momentum vector at the outermost hit position
Definition: Track.h:73
const TrackExtraRef & extra() const
reference to &quot;extra&quot; object
Definition: Track.h:97
HcalCalibrations calibped
const int ncidmx
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
Geom::Theta< T > theta() const
ReturnType plane(const PositionType &pos, const RotationType &rot) const
Definition: PlaneBuilder.h:22
std::vector< T >::const_iterator const_iterator
edm::ESHandle< HcalDbService > conditions_
T y() const
Definition: PV3DBase.h:63
#define abs(x)
Definition: mlp_lapack.h:159
std::map< std::string, bool > fired
#define min(a, b)
Definition: mlp_lapack.h:161
ErrorD< N >::type type
Definition: Error.h:30
double pedestal(int fCapId) const
get pedestal for capid=0..3
Geom::Phi< T > phi() const
const HcalQIEShape * m_shape
virtual void produce(edm::Event &, const edm::EventSetup &)
GlobalVector momentum() const
void findHOEtaPhi(int iphsect, int &ietaho, int &iphiho)
double charge(const std::vector< uint8_t > &Ampls)
Definition: Plane.h:17
void beginJob()
Definition: Breakpoints.cc:15
double outerZ() const
z coordinate of the outermost hit position
Definition: Track.h:81
static const unsigned int nL1mx
const math::XYZPoint & innerPosition() const
position of the innermost hit
Definition: Track.h:42
const int netahbmx
const HcalQIECoder * m_coder
int iEvent
Definition: GenABIO.cc:243
const T & max(const T &a, const T &b)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:94
float energy() const
Definition: CaloRecHit.h:19
GlobalPoint position() const
T sqrt(T t)
Definition: SSEVec.h:48
static const unsigned int nL1trg
Basic3DVector< float > RotationType
Vector3DBase< typename PreciseFloatType< T, U >::Type, FrameTag > cross(const Vector3DBase< U, FrameTag > &v) const
Definition: Vector3DBase.h:119
T z() const
Definition: PV3DBase.h:64
const int nsigpk
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
int ieta() const
get the cell ieta
Definition: HcalDetId.h:38
math::Error< 5 >::type CovarianceMatrix
int j
Definition: DBlmapReader.cc:9
Basic3DVector< float > PositionType
const int nstrbn
double outerX() const
x coordinate of the outermost hit position
Definition: Track.h:77
const int ntrgp_gm
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:361
float pedestal[netamx][nphimx][ncidmx]
double outerPz() const
z coordinate of momentum vector at the outermost hit position
Definition: Track.h:75
int k[5][pyjets_maxn]
int iphi() const
get the cell iphi
Definition: HcalDetId.h:40
TH1F * Nhotime[ntrgp_gm+1]
const T & get() const
Definition: EventSetup.h:55
std::vector< HOCalibVariables > HOCalibVariableCollection
collection of HOcalibration variabale
std::string const & label() const
Definition: InputTag.h:42
edm::EventID id() const
Definition: EventBase.h:56
AlCaHOCalibProducer(const edm::ParameterSet &)
T dot(const Basic3DVector &v) const
Scalar product, or &quot;dot&quot; product, with a vector of same type.
const int nphimx
const math::XYZVector & innerMomentum() const
momentum vector at the innermost hit position
Definition: Track.h:45
const int netamx
double outerY() const
y coordinate of the outermost hit position
Definition: Track.h:79
T * make() const
make new ROOT object
HcalCalibrationWidths calibwidth
int charge() const
track electric charge
Definition: TrackBase.h:113
const int nchnmx
TH1F * hotime[ntrgp_gm+1]
dbl *** dir
Definition: mlp_gen.cc:35
T x() const
Definition: PV3DBase.h:62
const int netahb3mx
static const unsigned int nHLTmx
HcalDetId id() const
get the id
Definition: HORecHit.h:21
double outerPx() const
x coordinate of momentum vector at the outermost hit position
Definition: Track.h:71
math::Error< dimension >::type CovarianceMatrix
5 parameter covariance matrix
Definition: TrackBase.h:72
T angle(T x1, T y1, T z1, T x2, T y2, T z2)
Definition: angle.h:11
float charge(const HcalQIEShape &fShape, unsigned fAdc, unsigned fCapId) const
ADC [0..127] + capid [0..3] -&gt; fC conversion.
Definition: HcalQIECoder.cc:22