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