CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/Calibration/HcalAlCaRecoProducers/src/AlCaHOCalibProducer.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Feb09 2009
00004 // Move the initialisation of SteppingHelixPropagator from ::beginJob() to ::produce()
00005 //
00006 // Oct3 2008
00007 // Difference in tag V00-02-45 with previous code
00008 
00009 // 1. One new object on data format, which was realised in  
00010 //     CRUZET data analysis.
00011 //2.  Remove all histogram and cout in the code
00012 //3. An upgrade in code, which increases the acceptance of 
00013 //    muon near the edge (this also realised in CRUZET data).
00014 // Difference in wrt V00-02-45
00015 // 1. initialisation tmpHOCalib.htime = -1000;
00016 // 2. By mistake HLT was commented out
00017 
00018 // Package:    AlCaHOCalibProducer
00019 // Class:      AlCaHOCalibProducer
00020 // 
00053 //
00054 // Original Author:  Gobinda Majumder
00055 //         Created:  Fri Jul  6 17:17:21 CEST 2007
00056 // $Id: AlCaHOCalibProducer.cc,v 1.27 2012/12/26 15:36:24 innocent Exp $
00057 //
00058 //
00059 
00060 
00061 // system include files
00062 #include <memory>
00063 
00064 // user include files
00065 #include "FWCore/Framework/interface/Frameworkfwd.h"
00066 #include "FWCore/Framework/interface/EDProducer.h"
00067 
00068 #include "FWCore/Framework/interface/Event.h"
00069 #include "FWCore/Framework/interface/MakerMacros.h"
00070 #include "FWCore/Framework/interface/EventSetup.h"
00071 #include "FWCore/Framework/interface/ESHandle.h"
00072 #include "DataFormats/FWLite/interface/Handle.h"
00073 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00074 
00075 #include "DataFormats/HcalDigi/interface/HcalDigiCollections.h"
00076 #include "FWCore/Utilities/interface/InputTag.h"
00077 #include "DataFormats/TrackReco/interface/Track.h"
00078 #include "DataFormats/TrackReco/interface/TrackFwd.h"
00079 
00080 #include "Geometry/Records/interface/IdealGeometryRecord.h"
00081 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00082 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00083 #include "Geometry/CaloGeometry/interface/CaloCellGeometry.h"
00084 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
00085 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
00086 
00087 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00088 #include "Geometry/DTGeometry/interface/DTGeometry.h"
00089 #include "Geometry/DTGeometry/interface/DTLayer.h"
00090 #include "Geometry/DTGeometry/interface/DTSuperLayer.h"
00091 #include "DataFormats/DTRecHit/interface/DTRecHitCollection.h"
00092 
00093 
00094 #include "CalibFormats/HcalObjects/interface/HcalDbService.h"
00095 #include "CalibFormats/HcalObjects/interface/HcalDbRecord.h"
00096 
00097 #include "CalibCalorimetry/HcalAlgos/interface/HcalAlgoUtils.h"
00098 #include "CalibCalorimetry/HcalAlgos/interface/HcalDbASCIIIO.h"
00099 #include "CalibFormats/HcalObjects/interface/HcalCalibrations.h"
00100 #include "CalibFormats/HcalObjects/interface/HcalCalibrationWidths.h"
00101 
00102 //08/07/07 #include "CondTools/Hcal/interface/HcalDbPool.h"
00103 //#include "CondFormats/HcalObjects/interface/HcalPedestals.h"
00104 //#include "CondFormats/HcalObjects/interface/HcalPedestalWidths.h"
00105 
00106 
00107 // #include "TrackingTools/GeomPropagators/interface/HelixArbitraryPlaneCrossing.h"
00108 #include "DataFormats/TrajectorySeed/interface/PropagationDirection.h"
00109 #include "DataFormats/GeometrySurface/interface/PlaneBuilder.h"
00110 
00111 #include "MagneticField/Engine/interface/MagneticField.h"
00112 #include "MagneticField/Records/interface/IdealMagneticFieldRecord.h"
00113 
00114 #include "DataFormats/HcalCalibObjects/interface/HOCalibVariables.h"
00115 #include "DataFormats/Math/interface/Error.h"
00116 #include "CLHEP/Vector/LorentzVector.h"
00117 
00118 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
00119 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
00120 
00121 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
00122 #include "DataFormats/Math/interface/Error.h"
00123 #include "TrackingTools/TrajectoryState/interface/FreeTrajectoryState.h"
00124 #include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixPropagator.h"
00125 //#include "TrackPropagation/SteppingHelixPropagator/interface/SteppingHelixStateInfo.h"
00126 
00127 #include "FWCore/ServiceRegistry/interface/Service.h"
00128 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00129 
00130 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00131 #include "DataFormats/Common/interface/TriggerResults.h"
00132 #include "FWCore/Common/interface/TriggerNames.h"
00133 
00134 
00135 #include "TFile.h"
00136 #include "TH1F.h"
00137 #include "TH2F.h"
00138 #include "TProfile.h"
00139 #include "TTree.h"
00140 /* C++ Headers */
00141 #include <string>
00142 
00143 #include <iostream>
00144 #include <fstream>
00145 //
00146 // class decleration
00147 //
00148 using namespace std;
00149 using namespace edm;
00150 using namespace reco;
00151 
00152 const int netabin= 16;  
00153 const int nphimx = 72;
00154 const int netamx = 32;
00155 const int nchnmx = 10;
00156 const int ncidmx = 5;
00157 
00158 //GMA #ifdef DEBUG_OFFLINE_GM
00159 const int nsigpk = 7;
00160 const int nstrbn = 0;
00161 const int ntrgp_gm = 11;
00162 
00163 
00164 const int netahbmx = 60;
00165 const int netahb3mx = 32;
00166 
00167 static const unsigned int nL1trg = 200;
00168 
00169 static const unsigned int nL1mx=140;
00170 static const unsigned int nHLTmx=140;
00171 //GMA #endif
00172 
00173 class AlCaHOCalibProducer : public edm::EDProducer {
00174    public:
00175       explicit AlCaHOCalibProducer(const edm::ParameterSet&);
00176       ~AlCaHOCalibProducer();
00177 
00178     typedef Basic3DVector<float>   PositionType;
00179     typedef Basic3DVector<float>   DirectionType;
00180     typedef Basic3DVector<float>   RotationType;
00181 
00182 
00183    private:
00184       void findHOEtaPhi(int iphsect, int& ietaho, int& iphiho);
00185       virtual void beginJob() ;
00186       virtual void produce(edm::Event&, const edm::EventSetup&);
00187       virtual void endJob() ;
00188 
00189       // ----------member data ---------------------------
00190 
00191   float  xhor0; //x-position in ring 0
00192   float  yhor0; //y-position in ring 0
00193   float  xhor1; //x-position in ring 1
00194   float  yhor1; //y-position in ring 1   
00195   int iring;    //Ring number -2,-1,0,1,2
00196 
00197   float  localxhor0; //local x-distance from edege in ring 0
00198   float  localyhor0; //local y-distance from edege in ring 0
00199   float  localxhor1; //local x-distance from edege in ring 1
00200   float  localyhor1; //local y-distance from edege in ring 1
00201 
00202   float pedestal[netamx][nphimx][ncidmx]; 
00203 
00204   std::string digiLabel;
00205   
00206   bool debug;
00207   std::string theRootFileName;
00208 
00209   //GMA #ifdef DEBUG_OFFLINE_GM
00210 
00211   TH1F* libhoped;
00212   TH1F* libhoped1;
00213 
00214   TH1F* allhotime;
00215   TH1F* hotime[ntrgp_gm+1];
00216   TH1F* hopedtime;
00217 
00218   TProfile* hopedpr;  
00219   TH1F* hopedrms;  
00220   TH1F* hst_hopedrms;    
00221 
00222   TProfile* hopeak[ntrgp_gm+1];
00223   TProfile* horatio;
00224 
00225   TH1F* Nallhotime;
00226   TH1F* Nhotime[ntrgp_gm+1];
00227   TH1F* Nhopedtime;
00228 
00229   TH1F* allhb1;
00230   TH1F* allhb2;
00231   TH1F* allhb3;
00232 
00233   TH1F* Nallhb1;
00234   TH1F* Nallhb2;
00235   TH1F* Nallhb3;
00236 
00237   TProfile* hb1pedpr;  
00238   TH1F* hb1pedrms;  
00239   TH1F* hst_hb1pedrms;    
00240 
00241   TH1F* ho_occupency[5];  
00242 
00243   bool m_hotime;
00244   //GM #endif
00245 
00246   edm::InputTag muonTags_;   // cosmicMuons or standAloneMuons
00247   edm::InputTag hbheLabel_;
00248   edm::InputTag hoLabel_;
00249   edm::InputTag hltLabel_;
00250   edm::InputTag l1Label_;  
00251   edm::InputTag towerLabel_;    
00252 
00253   bool m_digiInput;            // digi (true) or rechit (false)
00254   bool m_hbinfo;
00255   int m_startTS;
00256   int m_endTS;    
00257   double m_magscale;
00258   double m_sigma;
00259 
00260   typedef math::Error<5>::type CovarianceMatrix;  
00261   //#ifdef DEBUG_OFFLINE_GM
00262   //  int Nevents;
00263   int Noccu;
00264   //  int Npass;
00265   int nRuns;
00266   //#endif
00267 
00268   int irunold;
00269   //  SteppingHelixPropagator* stepProp;
00270   FreeTrajectoryState getFreeTrajectoryState( const reco::Track& tk, const MagneticField* field, int itag, bool dir);
00271 
00272   edm::ESHandle<HcalDbService> conditions_;
00273   const HcalQIEShape* m_shape;
00274   const HcalQIECoder* m_coder;
00275 
00276   HcalCalibrations calibped;
00277   HcalCalibrationWidths calibwidth;
00278 
00279   unsigned int Ntp; // # of HLT trigger paths (should be the same for all events!)
00280   std::map<std::string, bool> fired; 
00281 
00282 };
00283 
00284 //
00285 // constants, enums and typedefs
00286 //
00287 
00288 //
00289 // static data member definitions
00290 //
00291 
00292 //
00293 // constructors and destructor
00294 //
00295 AlCaHOCalibProducer::AlCaHOCalibProducer(const edm::ParameterSet& iConfig)
00296   :  muonTags_(iConfig.getUntrackedParameter<edm::InputTag>("muons"))
00297 
00298 {
00299    //register your products
00300 
00301   theRootFileName = iConfig.getUntrackedParameter<string>("RootFileName","tmp.root");
00302   m_digiInput = iConfig.getUntrackedParameter<bool>("digiInput", true);
00303   m_hbinfo = iConfig.getUntrackedParameter<bool>("hbinfo", false);
00304   m_startTS = iConfig.getUntrackedParameter<int>("firstTS", 4);
00305 
00306   m_hotime = iConfig.getUntrackedParameter<bool>("hotime", false);
00307 
00308   if(m_startTS<0) m_startTS=0;
00309   m_endTS = iConfig.getUntrackedParameter<int>("lastTS", 7);
00310   if (m_endTS < m_startTS) m_endTS = m_startTS + 3;
00311   if (m_endTS >9) m_endTS=9;
00312   m_magscale = iConfig.getUntrackedParameter<double>("m_scale", 4.0);
00313   m_sigma = iConfig.getUntrackedParameter<double>("sigma", 1.0);
00314   
00315   hoLabel_ = iConfig.getParameter<edm::InputTag>("hoInput");
00316   hbheLabel_ = iConfig.getParameter<edm::InputTag>("hbheInput");
00317   hltLabel_ = iConfig.getParameter<edm::InputTag>("hltInput");
00318   l1Label_ = iConfig.getParameter<edm::InputTag>("l1Input");
00319   towerLabel_ = iConfig.getParameter<edm::InputTag>("towerInput");  
00320   
00321   produces<HOCalibVariableCollection>("HOCalibVariableCollection").setBranchAlias("HOCalibVariableCollection");
00322   
00323   
00324   if (m_hotime) {
00325     edm::Service<TFileService> fs;
00326     
00327     char title[200];
00328     if ( m_digiInput) {
00329       libhoped = fs->make<TH1F>("libhoped", "libhoped", ncidmx*netamx*nphimx, -0.5, ncidmx*netamx*nphimx-0.5);
00330       libhoped1 = fs->make<TH1F>("libhoped1", "libhoped1", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00331       allhotime = fs->make<TH1F>("allhotime", "allhotime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00332       for (int ij=0; ij<=ntrgp_gm; ij++) {
00333         sprintf(title, "hotime_trgp_%i", ij+1);
00334         hotime[ij] = fs->make<TH1F>(title, title, nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00335         sprintf(title, "hopeak_trgp_%i", ij+1);
00336         hopeak[ij] = fs->make<TProfile>(title, title,netamx*nphimx, -0.5, netamx*nphimx-0.5);    
00337       }
00338       
00339       horatio = fs->make<TProfile>("horatio", "horatio",netamx*nphimx, -0.5, netamx*nphimx-0.5);    
00340       hopedtime = fs->make<TH1F>("hopedtime", "hopedtime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00341       
00342       Nallhotime = fs->make<TH1F>("Nallhotime", "Nallhotime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00343       hopedpr = fs->make<TProfile>("hopedpr", "hopedpr", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00344       hopedrms = fs->make<TH1F>("hopedrms", "hopedrms", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00345       hst_hopedrms = fs->make<TH1F>("hst_hopedrms", "hst_hopedrms", 100, 0.0, 0.1);
00346       for (int ij=0; ij<=ntrgp_gm; ij++) {
00347         sprintf(title, "Nhotime_trgp_%i", ij+1);
00348         Nhotime[ij] = fs->make<TH1F>(title, title, nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00349       }
00350       Nhopedtime = fs->make<TH1F>("Nhopedtime", "Nhopedtime", nchnmx*netamx*nphimx, -0.5, nchnmx*netamx*nphimx-0.5);
00351       allhb1 = fs->make<TH1F>("allhb1", "allhb1", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
00352       allhb2 = fs->make<TH1F>("allhb2", "allhb2", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5); 
00353       allhb3 = fs->make<TH1F>("allhb3", "allhb3", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5); 
00354       Nallhb1 = fs->make<TH1F>("Nallhb1", "Nallhb1", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
00355       Nallhb2 = fs->make<TH1F>("Nallhb2", "Nallhb2", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);
00356       Nallhb3 = fs->make<TH1F>("Nallhb3", "Nallhb3", nchnmx*netahb3mx*nphimx, -0.5, nchnmx*netahb3mx*nphimx-0.5);  
00357       hb1pedpr = fs->make<TProfile>("hb1pedpr", "hb1pedpr", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
00358       hb1pedrms = fs->make<TH1F>("hb1pedrms", "hb1pedrms", nchnmx*netahbmx*nphimx, -0.5, nchnmx*netahbmx*nphimx-0.5);
00359       hst_hb1pedrms = fs->make<TH1F>("hst_hb1pedrms", "hst_hb1pedrms", 100, 0., 0.1);
00360       
00361     }
00362     for (int i=0; i<5; i++) {
00363       sprintf(title, "ho_occupency (>%i #sigma)", i+2); 
00364       ho_occupency[i] = fs->make<TH1F>(title, title, netamx*nphimx, -0.5, netamx*nphimx-0.5); 
00365     }
00366   }
00367 
00368 }
00369 
00370 AlCaHOCalibProducer::~AlCaHOCalibProducer()
00371 {
00372  
00373   // do anything here that needs to be done at desctruction time
00374   // (e.g. close files, deallocate resources etc.)
00375 
00376   if (m_hotime) {
00377     //  Write the histos to file
00378     if ( m_digiInput) {
00379       allhotime->Divide(Nallhotime);
00380       for (int ij=0; ij<=ntrgp_gm; ij++) {
00381         hotime[ij]->Divide(Nhotime[ij]);
00382       }
00383       hopedtime->Divide(Nhopedtime);
00384       libhoped->Scale(1./max(1,nRuns));
00385       libhoped1->Scale(1./max(1,nRuns));   
00386       for (int i=0; i<nchnmx*netamx*nphimx; i++) {
00387         float xx = hopedpr->GetBinError(i+1);
00388         if (hopedpr->GetBinEntries(i+1) >0) {
00389           hopedrms->Fill(i, xx);
00390           hst_hopedrms->Fill(xx);
00391         }
00392       }
00393       allhb1->Divide(Nallhb1);
00394       allhb2->Divide(Nallhb2);
00395       allhb3->Divide(Nallhb3);
00396       for (int i=0; i<nchnmx*netahbmx*nphimx; i++) {
00397         float xx = hb1pedpr->GetBinError(i+1);
00398         if (hb1pedpr->GetBinEntries(i+1) >0) {
00399           hb1pedrms->Fill(i, xx);
00400           hst_hb1pedrms->Fill(xx);
00401         }
00402       }  
00403     }
00404     for (int i=0; i<5; i++) {
00405       ho_occupency[i]->Scale(1./max(1,Noccu));
00406     }
00407   }
00408 
00409 }
00410 
00411 
00412 //
00413 // member functions
00414 //
00415 
00416 // ------------ method called to produce the data  ------------
00417 void
00418 AlCaHOCalibProducer::produce(edm::Event& iEvent, const edm::EventSetup& iSetup)
00419 {
00420 
00421   using namespace edm;
00422   int irun = iEvent.id().run();
00423   if (m_digiInput) {
00424     if (irunold !=irun)  { 
00425       iSetup.get<HcalDbRecord>().get(conditions_);
00426 
00427       for (int i=0; i<netamx; i++) {
00428         for (int j=0; j<nphimx; j++) {
00429           for (int k=0; k<ncidmx; k++) {
00430             pedestal[i][j][k]=0.0;
00431           }
00432         }
00433       }     
00434     }
00435   }
00436 
00437   //  if (m_hotime && m_digiInput) {
00438   if (m_digiInput) {
00439     if (irunold !=irun) {
00440       nRuns++;
00441       for (int i =-netabin+1; i<=netabin-1; i++) {
00442         if (i==0) continue;
00443         int tmpeta1 =  (i>0) ? i -1 : -i +14; 
00444         if (tmpeta1 <0 || tmpeta1 >netamx) continue;
00445         for (int j=0; j<nphimx; j++) {
00446           
00447           HcalDetId id(HcalOuter, i, j+1, 4);
00448           calibped = conditions_->getHcalCalibrations(id);
00449           
00450           for (int k =0; k<ncidmx-1; k++) {
00451             pedestal[tmpeta1][j][k] = calibped.pedestal(k); // pedm->getValue(k);
00452             pedestal[tmpeta1][j][ncidmx-1] += (1./(ncidmx-1))*pedestal[tmpeta1][j][k];
00453           }
00454           
00455           if (m_hotime) {
00456             for (int k =0; k<ncidmx; k++) {
00457               libhoped->Fill(nphimx*ncidmx*tmpeta1 + ncidmx*j + k, pedestal[tmpeta1][j][k]);
00458             }
00459             for (int k =0; k<nchnmx; k++) {
00460               libhoped1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*j + k, pedestal[tmpeta1][j][min(k,ncidmx-1)]);
00461             }
00462           }
00463 
00464         }
00465       }
00466     }
00467   }
00468 
00469   //  Nevents++;
00470   irunold = irun;
00471 
00472   //GMA  if (Nevents%500==1) 
00473   //GMA  cout <<"AlCaHOCalibProducer Processing event # "<<Nevents<<" "<<Npass<<" "<<Noccu<<" "<<irun<<" "<<iEvent.id().event()<<endl;
00474 
00475   std::auto_ptr<HOCalibVariableCollection> hostore (new HOCalibVariableCollection);
00476 
00477   edm::Handle<HODigiCollection> ho;   
00478   
00479   edm::Handle<HBHEDigiCollection> hbhe; 
00480 
00481   if (m_digiInput) {
00482       iEvent.getByLabel(hoLabel_,ho);
00483       iEvent.getByLabel(hbheLabel_,hbhe);
00484   }
00485   
00486   if (m_hotime && m_digiInput) {
00487     if ((*ho).size()>0) {
00488       for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
00489         HcalDetId id =(*j).id();
00490         m_coder = (*conditions_).getHcalCoder(id);
00491         m_shape = (*conditions_).getHcalShape(m_coder);
00492         int tmpeta= id.ieta();
00493         int tmpphi= id.iphi();
00494         float tmpdata[nchnmx];
00495         int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
00496         for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00497           tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00498           allhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
00499           Nallhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.);
00500         }
00501       }
00502     }
00503     if ((*hbhe).size()>0) {
00504       for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
00505         HcalDetId id =(*j).id();
00506         m_coder = (*conditions_).getHcalCoder(id);
00507         m_shape = (*conditions_).getHcalShape(m_coder);
00508         int tmpeta= id.ieta();
00509         int tmpphi= id.iphi();
00510         int tmpdepth =id.depth();
00511         int tmpeta1 =  (tmpeta>0) ? tmpeta -15 : -tmpeta + 1; 
00512         if (tmpdepth==1) tmpeta1 =  (tmpeta>0) ? tmpeta -1 : -tmpeta +29;  
00513         for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00514           float signal = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00515           if (tmpdepth==1) { 
00516             allhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00517             Nallhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);
00518             hb1pedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);}
00519           if (tmpdepth==2) { 
00520             allhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00521             Nallhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
00522           if (tmpdepth==3) { 
00523             allhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00524             Nallhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
00525         }
00526       }
00527     }
00528   }
00529 
00530   double pival = acos(-1.);
00531   
00532   Handle<reco::TrackCollection> cosmicmuon;
00533   iEvent.getByLabel(muonTags_, cosmicmuon);
00534   
00535   if (cosmicmuon->size()>0) { 
00536     
00537     int l1trg = 0;
00538     int hlttr = 0;
00539     
00540     int ntrgpas_gm[ntrgp_gm]={0,0,0,0,0,0,0,0,0,0};
00541  
00542     /*   
00543     //L1 trigger
00544     Handle<L1GlobalTriggerReadoutRecord> L1GTRR;
00545     iEvent.getByLabel(l1Label_,L1GTRR);  //gtDigis
00546     
00547     if ( L1GTRR.isValid()) {
00548       const unsigned int n(L1GTRR->decisionWord().size());
00549       const bool accept(L1GTRR->decision());
00550       if (accept) {
00551         for (unsigned int i=0; i!=n && i<32; ++i) {
00552           //    for (unsigned int i=0; i!=n ; ++i) {
00553           int il1trg = (L1GTRR->decisionWord()[i]) ? 1 : 0;
00554           if (il1trg>0 && i<32) l1trg +=int(std::pow(2., double(i%32))*il1trg);
00555         }
00556       }
00557     }// else { return;}
00558     
00559     //HLT 
00560 
00561     Handle<edm::TriggerResults> trigRes;    
00562     iEvent.getByLabel(hltLabel_, trigRes);
00563 
00564 
00565     unsigned int size = trigRes->size();
00566     edm::TriggerNames triggerNames(*trigRes);
00567     
00568     // loop over all paths, get trigger decision
00569     for(unsigned i = 0; i != size && i<32; ++i) {
00570       std::string name = triggerNames.triggerName(i);
00571       fired[name] = trigRes->accept(i);
00572       int ihlt =  trigRes->accept(i);
00573       if (m_hotime){ 
00574         if (ihlt >0 && i < (int)ntrgp_gm) { ntrgpas_gm[i] = 1;}
00575       }
00576       if (i<32 && ihlt>0) hlttr += int(std::pow(2., double(i%32))*ihlt);
00577     }
00578 
00579     */
00580 
00581     int Noccu_old = Noccu;
00582     
00583     for(reco::TrackCollection::const_iterator ncosm = cosmicmuon->begin();
00584         ncosm != cosmicmuon->end();  ++ncosm) {
00585       
00586       if ((*ncosm).ndof() < 15) continue;
00587       if ((*ncosm).normalizedChi2() >30.0) continue;
00588 
00589       HOCalibVariables tmpHOCalib;
00590       
00591       tmpHOCalib.trig1 = l1trg;
00592       tmpHOCalib.trig2 = hlttr;    
00593       
00594       int charge = ncosm->charge();  
00595       
00596       double innerr = (*ncosm).innerPosition().Perp2();
00597       double outerr = (*ncosm).outerPosition().Perp2();
00598       int iiner = (innerr <outerr) ? 1 : 0;
00599       
00600       //---------------------------------------------------
00601       //             in_to_out  Dir         in_to_out  Dir
00602       //   StandAlone ^         ^     Cosmic    ^    |
00603       //              |         |               |    v
00604       //---------------------------------------------------Y=0
00605       //   StandAlone |         |     Cosmic    ^    |
00606       //              v         v               |    v
00607       //----------------------------------------------------
00608       
00609       double posx, posy, posz;
00610       double momx, momy, momz;
00611       
00612       if (iiner==1) {
00613         posx = (*ncosm).innerPosition().X();
00614         posy = (*ncosm).innerPosition().Y();
00615         posz = (*ncosm).innerPosition().Z();
00616         
00617         momx = (*ncosm).innerMomentum().X();
00618         momy = (*ncosm).innerMomentum().Y();
00619         momz = (*ncosm).innerMomentum().Z();
00620         
00621       } else {
00622         posx = (*ncosm).outerPosition().X();
00623         posy = (*ncosm).outerPosition().Y();
00624         posz = (*ncosm).outerPosition().Z();
00625         
00626         momx = (*ncosm).outerMomentum().X();
00627         momy = (*ncosm).outerMomentum().Y();
00628         momz = (*ncosm).outerMomentum().Z();
00629       }
00630       
00631       
00632       PositionType trkpos(posx, posy, posz);
00633       
00634       CLHEP::Hep3Vector tmpmuon3v(posx, posy, posz);
00635       CLHEP::Hep3Vector tmpmuondir(momx, momy, momz);
00636       
00637       bool samedir = (tmpmuon3v.dot(tmpmuondir) >0) ? true : false;
00638       for (int i=0; i<3; i++) {tmpHOCalib.caloen[i] = 0.0;}
00639       int inearbymuon = 0;
00640       for(reco::TrackCollection::const_iterator ncosmcor = cosmicmuon->begin();
00641           ncosmcor != cosmicmuon->end();  ++ncosmcor) {
00642         if (ncosmcor==ncosm) continue;
00643         
00644         CLHEP::Hep3Vector tmpmuon3vcor;
00645         CLHEP::Hep3Vector tmpmom3v;
00646         if (iiner==1) {
00647           tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).innerPosition().X(),(*ncosmcor).innerPosition().Y(),(*ncosmcor).innerPosition().Z());
00648           tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).innerMomentum().X(),(*ncosmcor).innerMomentum().Y(),(*ncosmcor).innerMomentum().Z());
00649         } else {
00650           tmpmuon3vcor = CLHEP::Hep3Vector((*ncosmcor).outerPosition().X(),(*ncosmcor).outerPosition().Y(),(*ncosmcor).outerPosition().Z());
00651           tmpmom3v = CLHEP::Hep3Vector((*ncosmcor).outerMomentum().X(),(*ncosmcor).outerMomentum().Y(),(*ncosmcor).outerMomentum().Z());        
00652           
00653         }
00654         if (tmpmom3v.mag()<0.2 || (*ncosmcor).ndof()<5) continue;
00655         
00656         double angle = tmpmuon3v.angle(tmpmuon3vcor);
00657         if (angle < 7.5*pival/180.) {inearbymuon=1;} //  break;}
00658 
00659         if (muonTags_.label() =="cosmicMuons") {
00660           if (angle <7.5*pival/180.) { tmpHOCalib.caloen[0] +=1.;}
00661           if (angle <15.0*pival/180.) { tmpHOCalib.caloen[1] +=1.;}
00662           if (angle <35.0*pival/180.) { tmpHOCalib.caloen[2] +=1.;}
00663         }
00664       }
00665       
00666       localxhor0 = localyhor0 = 20000;  //GM for 22OCT07 data
00667       
00668       if (muonTags_.label() =="standAloneMuons") {
00669         
00670         Handle<CaloTowerCollection> calotower;
00671         iEvent.getByLabel(towerLabel_, calotower);
00672 
00673         for (CaloTowerCollection::const_iterator calt = calotower->begin();
00674              calt !=calotower->end(); calt++) {
00675           //CMSSW_2_1_x const math::XYZVector towermom = (*calt).momentum();
00676           double ith = (*calt).momentum().theta();
00677           double iph = (*calt).momentum().phi();
00678           
00679           CLHEP::Hep3Vector calo3v(sin(ith)*cos(iph), sin(ith)*sin(iph), cos(ith));
00680           
00681           double angle = tmpmuon3v.angle(calo3v);
00682           
00683           if (angle < 7.5*pival/180.) {tmpHOCalib.caloen[0] += calt->emEnergy()+calt->hadEnergy();}
00684           if (angle < 15*pival/180.) {tmpHOCalib.caloen[1] += calt->emEnergy()+calt->hadEnergy();}
00685           if (angle < 35*pival/180.) {tmpHOCalib.caloen[2] += calt->emEnergy()+calt->hadEnergy();}
00686         }
00687         
00688         
00689       }
00690       if (tmpHOCalib.caloen[0] >10.0) continue;
00691       
00692       GlobalPoint glbpt(posx, posy, posz);
00693       
00694       double mom = sqrt(momx*momx + momy*momy +momz*momz);
00695       
00696       momx /= mom;
00697       momy /= mom;
00698       momz /= mom;
00699       
00700       DirectionType trkdir(momx, momy, momz);
00701       
00702       tmpHOCalib.trkdr = (*ncosm).d0();
00703       tmpHOCalib.trkdz = (*ncosm).dz();
00704       
00705       tmpHOCalib.nmuon = cosmicmuon->size();
00706       tmpHOCalib.trkvx = glbpt.x();
00707       tmpHOCalib.trkvy = glbpt.y();
00708       tmpHOCalib.trkvz = glbpt.z();
00709       tmpHOCalib.trkmm = mom*charge;
00710       tmpHOCalib.trkth = trkdir.theta();
00711       tmpHOCalib.trkph = trkdir.phi();
00712       
00713       tmpHOCalib.ndof  = (inearbymuon ==0) ? (int)(*ncosm).ndof() : -(int)(*ncosm).ndof();
00714       tmpHOCalib.chisq = (*ncosm).normalizedChi2(); // max(1.,tmpHOCalib.ndof);
00715       tmpHOCalib.therr = 0.;
00716       tmpHOCalib.pherr = 0.;
00717       
00718       if (iiner==1) {
00719         reco::TrackBase::CovarianceMatrix innercov = (*ncosm).innerStateCovariance();
00720         tmpHOCalib.therr = innercov(1,1); //thetaError();
00721         tmpHOCalib.pherr = innercov(2,2); //phi0Error();
00722       } else {
00723         reco::TrackBase::CovarianceMatrix outercov = (*ncosm).outerStateCovariance();
00724         tmpHOCalib.therr = outercov(1,1); //thetaError();
00725         tmpHOCalib.pherr = outercov(2,2); //phi0Error();
00726       }
00727       
00728       ESHandle<MagneticField> theMagField;
00729       iSetup.get<IdealMagneticFieldRecord>().get(theMagField );     
00730 
00731       SteppingHelixPropagator myHelix(&*theMagField,anyDirection);
00732       myHelix.setMaterialMode(false);
00733       myHelix.applyRadX0Correction(true);
00734 
00735       double phiho = trkpos.phi();
00736       if (phiho<0) phiho +=2*pival;
00737       
00738       int iphisect_dt=int(6*(phiho+pival/18.)/pival); //for u 18/12/06
00739       if (iphisect_dt>=12) iphisect_dt=0;
00740 
00741       int iphisect = -1;
00742 
00743       bool ipath = false;
00744       for (int kl = 0; kl<=2; kl++) {
00745 
00746         int iphisecttmp = (kl<2) ? iphisect_dt + kl : iphisect_dt - 1;
00747         if (iphisecttmp <0) iphisecttmp = 11;
00748         if (iphisecttmp >=12) iphisecttmp = 0;
00749         
00750         double phipos = iphisecttmp*pival/6.;
00751         double phirot = phipos;
00752         
00753         GlobalVector xLocal(-sin(phirot), cos(phirot), 0.);
00754         
00755         GlobalVector yLocal(0., 0., 1.);
00756         GlobalVector zLocal = xLocal.cross(yLocal).unit();
00757         //    GlobalVector zLocal(cos(phirot), sin(phirot), 0.0); 
00758         
00759 
00760         FreeTrajectoryState freetrajectorystate_ = getFreeTrajectoryState(*ncosm,&(*theMagField), iiner, samedir);
00761         
00762         Surface::RotationType rot(xLocal, yLocal, zLocal);
00763         
00764         for (int ik=1; ik>=0; ik--) { //propagate track in two HO layers
00765           
00766           double radial = 407.0;
00767           if (ik==0) radial = 382.0;
00768 
00769           Surface::PositionType pos(radial*cos(phipos), radial*sin(phipos), 0.);
00770           PlaneBuilder::ReturnType aPlane = PlaneBuilder().plane(pos,rot);
00771 
00772           auto aPlane2 = new Plane(pos,rot);
00773 
00774           SteppingHelixStateInfo steppingHelixstateinfo_ = myHelix.propagate(SteppingHelixStateInfo(freetrajectorystate_), (*aPlane2));
00775 
00776           if (steppingHelixstateinfo_.isValid()) {
00777 
00778             GlobalVector hotrkpos2(steppingHelixstateinfo_.position().x(), steppingHelixstateinfo_.position().y(), steppingHelixstateinfo_.position().z());
00779             CLHEP::Hep3Vector hotrkdir2(steppingHelixstateinfo_.momentum().x(), steppingHelixstateinfo_.momentum().y(),steppingHelixstateinfo_.momentum().z());
00780             
00781             LocalVector lclvt0 = (*aPlane).toLocal(hotrkpos2);
00782             
00783             double xx = lclvt0.x();
00784             double yy = lclvt0.y();
00785             
00786             if (ik ==1) {
00787               if ((std::abs(yy) < 130 && xx >-64.7 && xx <138.2)
00788                   ||(std::abs(yy) > 130 && std::abs(yy) <700 && xx >-76.3 && xx <140.5)) {
00789                 ipath = true;  //Only look for tracks which as hits in layer 1
00790                 iphisect = iphisecttmp;
00791               }
00792             }
00793             
00794             if (iphisect != iphisecttmp) continue; //Look for ring-0 only when ring1 is accepted for that sector
00795             
00796             switch (ik) 
00797               {
00798               case 0 : 
00799                 xhor0 = xx; //lclvt0.x();
00800                 yhor0 = yy; //lclvt0.y();
00801                 break;
00802               case 1 :
00803                 xhor1 = xx; //lclvt0.x();
00804                 yhor1 = yy; //lclvt0.y();
00805                 
00806                 tmpHOCalib.hoang = CLHEP::Hep3Vector(zLocal.x(),zLocal.y(),zLocal.z()).dot(hotrkdir2.unit());
00807                 break;
00808               default : break;
00809               }
00810           } else {
00811             break;
00812           }
00813         }
00814         if (ipath) break;
00815       }
00816       if (ipath) { //If muon crossed HO laeyrs
00817         
00818         int ietaho = 50;
00819         int iphiho = -1;
00820         
00821         for (int i=0; i<9; i++) {tmpHOCalib.hosig[i]=-100.0;}
00822         for (int i=0; i<18; i++) {tmpHOCalib.hocorsig[i]=-100.0;}
00823         for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
00824         tmpHOCalib.hocro = -100;
00825         tmpHOCalib.htime = -1000;
00826         
00827         int isect = 0;
00828 
00829         findHOEtaPhi(iphisect, ietaho, iphiho);
00830         
00831         if (ietaho !=0 && iphiho !=0 && std::abs(iring)<=2) { //Muon passed through a tower
00832           isect = 100*std::abs(ietaho+30)+std::abs(iphiho);
00833           if (std::abs(ietaho) >=netabin || iphiho<0) isect *=-1; //Not extrapolated to any tower
00834           if (std::abs(ietaho) >=netabin) isect -=1000000;  //not matched with eta
00835           if (iphiho<0)        isect -=2000000; //not matched with phi
00836           tmpHOCalib.isect = isect;
00837           
00838           tmpHOCalib.hodx = localxhor1;
00839           tmpHOCalib.hody = localyhor1;      
00840           
00841           if (iring==0) {
00842             tmpHOCalib.hocorsig[8] = localxhor0;
00843             tmpHOCalib.hocorsig[9] = localyhor0;
00844           }
00845           
00846           int etamn=-4;
00847           int etamx=4;
00848           if (iring==1) {etamn=5; etamx = 10;}
00849           if (iring==2) {etamn=11; etamx = 16;}
00850           if (iring==-1){etamn=-10; etamx = -5;}
00851           if (iring==-2){etamn=-16; etamx = -11;}
00852           
00853           int phimn = 1;
00854           int phimx = 2;
00855           if (iring ==0) {
00856             phimx =2*int((iphiho+1)/2.);
00857             phimn = phimx - 1;
00858           } else {
00859             phimn = 3*int((iphiho+1)/3.) - 1; 
00860             phimx = phimn + 2;
00861           }
00862           
00863           if (phimn <1) phimn += nphimx;
00864           if (phimx >72) phimx -= nphimx;
00865           
00866           int sigstr = m_startTS; // 5;
00867           int sigend = m_endTS; // 8;
00868           
00869           //      if (iphiho <=nphimx/2) { //GMA310508
00870           //        sigstr -=1; //GMA300608
00871           //        sigend -=1;
00872           //      }
00873           
00874           if (m_hbinfo) {
00875             for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
00876             
00877             if (m_digiInput) {
00878               if ((*hbhe).size() >0) {
00879                 for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
00880                   //              const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00881                   //              HcalDetId id =digi.id();
00882                   HcalDetId id =(*j).id();
00883                   m_coder = (*conditions_).getHcalCoder(id);
00884                   m_shape = (*conditions_).getHcalShape(m_coder);
00885                   int tmpeta= id.ieta();
00886                   int tmpphi= id.iphi();
00887                   calibped = conditions_->getHcalCalibrations(id);
00888                   
00889                   int deta = tmpeta-ietaho;
00890                   if (tmpeta==-1 && ietaho== 1) deta = -1;
00891                   if (tmpeta== 1 && ietaho==-1) deta =  1;
00892                   int dphi = tmpphi-iphiho;
00893                   if (phimn >phimx) {
00894                     if (dphi==71) dphi=-1;
00895                     if (dphi==-71) dphi=1;
00896                   }
00897                   
00898                   int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
00899                   
00900                   if (ipass2 ==0 ) continue;
00901                   
00902                   float tmpdata[nchnmx];
00903                   for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00904                     tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00905                   }
00906                   
00907                   float signal = 0;
00908                   for (int i=1; i<(*j).size() && i<=8; i++) {
00909                     signal += tmpdata[i] - calibped.pedestal((*j).sample(i).capid());; 
00910                   }
00911                   
00912                   if (ipass2 == 1) {
00913                     if (3*(deta+1)+dphi+1<9)  tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
00914                   }
00915                 }
00916               }
00917               
00918             } else {
00919               
00920               edm::Handle<HBHERecHitCollection> hbheht;// iEvent.getByType(hbheht);
00921               iEvent.getByLabel(hbheLabel_,hbheht);
00922 
00923               
00924               if ((*hbheht).size()>0) {
00925                 if(!(*hbheht).size()) throw (int)(*hbheht).size();
00926                 
00927                 for (HBHERecHitCollection::const_iterator j=(*hbheht).begin(); j!=(*hbheht).end(); j++){
00928                   //              const HBHERecHit hbhehtrec = (const HBHERecHit)(*j);
00929                   //              HcalDetId id =hbhehtrec.id();
00930                   HcalDetId id =(*j).id();
00931                   int tmpeta= id.ieta();
00932                   int tmpphi= id.iphi();
00933                   
00934                   int deta = tmpeta-ietaho;
00935                   if (tmpeta==-1 && ietaho== 1) deta = -1;
00936                   if (tmpeta== 1 && ietaho==-1) deta =  1;
00937                   int dphi = tmpphi-iphiho;
00938                   if (phimn >phimx) {
00939                     if (dphi==71) dphi=-1;
00940                     if (dphi==-71) dphi=1;
00941                   }
00942                   
00943                   int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
00944                   if ( ipass2 ==0 ) continue;
00945                   
00946                   float signal = (*j).energy();
00947                   
00948                   if (3*(deta+1)+dphi+1<9)  tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
00949                 }
00950               }
00951               
00952             } //else m_digilevel
00953             
00954           } //m_hbinfo #endif
00955           
00956           if (m_digiInput) {
00957             if ((*ho).size()>0) {
00958               int isFilled[netamx*nphimx]; 
00959               for (int j=0; j<netamx*nphimx; j++) {isFilled[j]=0;}
00960               
00961               double sumEt = 0;
00962               double sumE  = 0;
00963               
00964               for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
00965                 //              const HODataFrame digi = (const HODataFrame)(*j);
00966                 //              HcalDetId id =digi.id();
00967 
00968                 HcalDetId id =(*j).id();                
00969                 m_coder = (*conditions_).getHcalCoder(id);
00970                 m_shape = (*conditions_).getHcalShape(m_coder);
00971                 int tmpeta= id.ieta();
00972                 int tmpphi= id.iphi();
00973                 
00974                 int ipass1 =0;
00975                 if (tmpeta >=etamn && tmpeta <=etamx) {
00976                   if (phimn < phimx) {
00977                     ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
00978                   } else {
00979                     ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
00980                   }
00981                 }
00982                 
00983                 int deta = tmpeta-ietaho;
00984                 if (tmpeta==-1 && ietaho== 1) deta = -1;
00985                 if (tmpeta== 1 && ietaho==-1) deta =  1;
00986                 
00987                 int dphi = tmpphi -iphiho;
00988                 if (phimn>phimx) {
00989                   if (dphi==71) dphi=-1;
00990                   if (dphi==-71) dphi=1;
00991                 }
00992                 
00993                 int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0;
00994                 
00995                 int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
00996                 
00997                 float tmpdata[nchnmx]={0,0,0,0,0,0,0,0,0,0};
00998                 float sigvall[nsigpk]={0,0,0,0,0,0,0};
00999              
01000                 for (int i=0; i<(*j).size() && i<nchnmx; i++) {
01001                   tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
01002                   if (deta==0 && dphi==0) { 
01003                     double tmpE = tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01004                     if (tmpE >0) {
01005                       sumEt += i*tmpE;
01006                       sumE  += tmpE;
01007                     }
01008                     if (m_hotime) {
01009                       //calculate signals in 4 time slices, 0-3,.. 6-9
01010                       if (i>=7-nsigpk) {
01011                         for (int ncap=0; ncap<nsigpk; ncap++) {
01012                           if (i-ncap >= nstrbn && i-ncap <= nstrbn+3) { 
01013                             sigvall[ncap] +=tmpdata[i];
01014                           }
01015                         }
01016                         }
01017                       if (i==(*j).size()-1) {
01018                         float mxled=-1;
01019                         int imxled = 0;
01020                         for (int ij=0; ij<nsigpk; ij++) {
01021                           if (sigvall[ij] > mxled) {mxled = sigvall[ij]; imxled=ij;}
01022                         }
01023                         double pedx = 0.0;
01024                         for (int ij=0; ij<4; ij++) {
01025                           pedx +=pedestal[tmpeta1][tmpphi-1][ij];
01026                         }
01027                         if (mxled-pedx >2 && mxled-pedx <20 ) {
01028                           hopeak[ntrgp_gm]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
01029                           for (int jk=0; jk<ntrgp_gm; jk++) {
01030                             if (ntrgpas_gm[jk]>0) {
01031                               hopeak[jk]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
01032                             }
01033                           }
01034                           if (tmpdata[5]+tmpdata[6] >1) {
01035                             horatio->Fill(nphimx*tmpeta1 + tmpphi-1, (tmpdata[5]-tmpdata[6])/(tmpdata[5]+tmpdata[6]));
01036                           }
01037                           for (int ij=0; ij<(*j).size() && ij<nchnmx; ij++) {
01038                             hotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
01039                             Nhotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
01040                             for (int jk=0; jk<ntrgp_gm; jk++) {
01041                               if (ntrgpas_gm[jk]>0) {
01042                                 hotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
01043                                 Nhotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
01044                               }
01045                             }
01046                           }
01047                         }
01048                       }
01049                     }
01050                   }
01051                 }
01052 
01053                 if (std::abs(tmpeta) <=15 && deta==0 && dphi ==0) { 
01054                   float signal = 0;
01055                   size_t icnt = 0;
01056                   for (int i =0; i<nchnmx && i< (*j).size(); i++) {
01057                     if (i >=sigstr && i<=sigend) continue;
01058                     signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01059                     if (++icnt >=4) break;
01060                   }
01061                   tmpHOCalib.hocro = signal;
01062                 }
01063                 
01064                 if (m_hotime) { 
01065                   if (ipass1 ==0 && ipass2 ==0 && cosmicmuon->size()<=2) {
01066                     if (std::abs(ietaho) <=netabin && iphiho >0) {
01067                       if ((iphiho >=1 && iphiho<=nphimx/2 && tmpphi >=1 && tmpphi <=nphimx/2) ||
01068                           (iphiho >nphimx/2 && iphiho<=nphimx && tmpphi >nphimx/2 && tmpphi <=nphimx)) {
01069                         if (isFilled[nphimx*tmpeta1+tmpphi-1]==0) {
01070                           isFilled[nphimx*tmpeta1+tmpphi-1]=1;
01071                           for (int i=0; i<(*j).size() && i<nchnmx; i++) {
01072                             hopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
01073                             Nhopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.); 
01074                             hopedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
01075                           }
01076                         } //isFilled
01077                       }
01078                     }
01079                   }
01080                 }
01081                 
01082                 if (ipass1 ==0 && ipass2 ==0 ) continue;
01083                 
01084                 float signal = 0;
01085                 for (int i=sigstr; i<(*j).size() && i<=sigend; i++) {
01086                   signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01087                 }
01088                 if (signal <-100 || signal >100000) signal = -100;
01089                 if (m_hotime) {
01090                   if (signal >-100 && Noccu == Noccu_old) {
01091                     for (int i=0; i<5; i++) {
01092                       if (signal >(i+2)*m_sigma) {
01093                         ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
01094                       }
01095                     }
01096                   }
01097                 }
01098 
01099                 if (ipass1 ==0 && ipass2 ==0 ) continue;
01100                 
01101                 if (ipass1 ==1) {
01102                   int tmpdph = tmpphi-phimn;
01103                   if (tmpdph<0) tmpdph = 2;  //only case of iphi==1, where phimn=71
01104                   
01105                   int ilog = 2*(tmpeta-etamn)+tmpdph;
01106                   if (iring !=0) { 
01107                     if (iring >0) {
01108                       ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
01109                     } else {
01110                       ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
01111                     }
01112                   }
01113                   if (ilog>-1 && ilog<18) { 
01114                     tmpHOCalib.hocorsig[ilog] = signal;
01115                   }
01116                 }             
01117                 
01118                 if (ipass2 ==1) {
01119                   if (3*(deta+1)+dphi+1<9) tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
01120                 }
01121                 
01122                 /*
01123                 // Possibility to store pedestal by shifting phi tower by 6
01124                 // But, due to missing tower at +-5, we do not have always proper
01125                 // statistics and also in pedestal subtracted data, we do not have
01126                 // signal in that tower
01127                 // 
01128                 if (deta==0 && dphi ==0) {
01129                   int crphi = tmpphi + 6;
01130                   if (crphi >72) crphi -=72;
01131                   
01132                   for (HODigiCollection::const_iterator jcr=(*ho).begin(); jcr!=(*ho).end(); jcr++){
01133                   //                const HODataFrame (*jcr) = (const HODataFrame)(*jcr);
01134                   //                HcalDetId idcr =(*jcr).id();
01135                   HcalDetId id =(*jcr).id();
01136                   m_coder = (*conditions_).getHcalCoder(idcr);
01137                   m_shape = (*conditions_).getHcalShape(m_coder);
01138                     int etacr= idcr.ieta();
01139                     int phicr= idcr.iphi();
01140                     
01141                     if (tmpeta==etacr && crphi ==phicr) {
01142                       
01143                       float tmpdatacr[nchnmx];
01144                       for (int i=0; i<(*jcr).size() && i<nchnmx; i++) {
01145                         tmpdatacr[i] = m_coder->charge(*m_shape,(*jcr).sample(i).adc(),(*jcr).sample(i).capid());
01146                       }
01147                     }
01148                     }
01149                     }
01150                 */
01151                 
01152             }
01153             tmpHOCalib.htime = sumEt/max(sumE,1.e-6);
01154           } 
01155         } else {
01156           edm::Handle<HORecHitCollection> hoht;
01157           iEvent.getByLabel(hoLabel_,hoht);
01158             
01159           
01160           if ((*hoht).size()>0) {
01161             for (HORecHitCollection::const_iterator j=(*hoht).begin(); j!=(*hoht).end(); j++){
01162               //                const HORecHit hohtrec = (const HORecHit)(*j);
01163               //                HcalDetId id =hohtrec.id();
01164               HcalDetId id =(*j).id();
01165               int tmpeta= id.ieta();
01166               int tmpphi= id.iphi();
01167 
01168               int ipass1 =0;
01169               if (tmpeta >=etamn && tmpeta <=etamx) {
01170                 if (phimn < phimx) {
01171                   ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
01172                 } else {
01173                   ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
01174                 }
01175               }
01176               
01177               int deta = tmpeta-ietaho;
01178               if (tmpeta==-1 && ietaho== 1) deta = -1;
01179               if (tmpeta== 1 && ietaho==-1) deta =  1;
01180               
01181               int dphi = tmpphi -iphiho;
01182               if (phimn>phimx) {
01183                 if (dphi==71) dphi=-1;
01184                 if (dphi==-71) dphi=1;
01185               }
01186               
01187               float signal = (*j).energy();
01188               if (m_hotime) {
01189                 int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
01190                 if (signal >-100 && Noccu == Noccu_old) {
01191                   for (int i=0; i<5; i++) {
01192                     if (signal >(i+2)*m_sigma) {
01193                       ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
01194                     }
01195                   }
01196                 }
01197               }
01198               
01199               int ipass2 = (std::abs(deta) <=1 && std::abs(dphi)<=1) ? 1 : 0;
01200               
01201               if (ipass1 ==0 && ipass2 ==0 ) continue;
01202               
01203               if (ipass1 ==1) {
01204                 int tmpdph = tmpphi-phimn;
01205                 if (tmpdph<0) tmpdph = 2;  //only case of iphi==1, where phimn=71
01206                   
01207                 int ilog = 2*(tmpeta-etamn)+tmpdph;
01208                 if (iring !=0) { 
01209                   if (iring >0) {
01210                     ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
01211                   } else {
01212                     ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
01213                   }
01214                 }
01215                 if (ilog>-1 && ilog<18) {
01216                   tmpHOCalib.hocorsig[ilog] = signal;
01217                 }
01218               }       
01219               
01220               if (ipass2 ==1) {
01221                 
01222                 if (3*(deta+1)+dphi+1<9) {
01223                   tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
01224                 }
01225               }
01226               
01227               if (deta==0 && dphi ==0) {
01228                 tmpHOCalib.htime = (*j).time();
01229                 int crphi = tmpphi + 6;
01230                 if (crphi >72) crphi -=72;
01231                 
01232                 for (HORecHitCollection::const_iterator jcr=(*hoht).begin(); jcr!=(*hoht).end(); jcr++){
01233                   const HORecHit reccr = (const HORecHit)(*jcr);
01234                   HcalDetId idcr =reccr.id();
01235                   int etacr= idcr.ieta();
01236                   int phicr= idcr.iphi();
01237                   if (tmpeta==etacr && crphi ==phicr) {
01238                     
01239                     tmpHOCalib.hocro = reccr.energy();
01240                     
01241                   }
01242                 }
01243               }
01244             }
01245           } 
01246         }
01247         
01248         //GMA     Npass++;
01249         if (Noccu == Noccu_old) Noccu++;
01250         hostore->push_back(tmpHOCalib); 
01251         
01252       }
01253     }
01254     
01255     } 
01256   } 
01257 
01258   iEvent.put(hostore, "HOCalibVariableCollection");
01259   
01260 }
01261 
01262 // ------------ method called once each job just before starting event loop  ------------
01263 void 
01264 AlCaHOCalibProducer::beginJob()
01265 {
01266   //GMA  Nevents = 0;
01267   //GMA  Npass = 0;
01268   //GMA  Noccu = 0;
01269 
01270   irunold = -1;
01271   nRuns = 0;
01272   //  edm::ESHandle<MagneticField> bField;
01273   //  iSetup.get<IdealMagneticFieldRecord>().get(bField);
01274   //  stepProp  = new SteppingHelixPropagator(&*bField,anyDirection);
01275   //  stepProp->setMaterialMode(false);
01276   //  stepProp->applyRadX0Correction(true);
01277   
01278   for (int i=0; i<netamx; i++) {
01279     for (int j=0; j<nphimx; j++) {
01280       for (int k=0; k<ncidmx; k++) {
01281         pedestal[i][j][k]=0.0;
01282       }
01283     }
01284   }
01285 
01286 
01287 }
01288 
01289 // ------------ method called once each job just after ending the event loop  ------------
01290 void 
01291 AlCaHOCalibProducer::endJob() {
01292 
01293 
01294 }
01295 
01296 void AlCaHOCalibProducer::findHOEtaPhi(int iphisect, int& ietaho, int& iphiho) {
01297   
01298   //18/12/06 : use only position, not angle phi
01299 
01300 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};
01301 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};
01302 
01303   double philow[6]={-76.27, -35.11, 0.35, 35.81, 71.77, 108.93};  //Ring+/-1 & 2
01304   double phihgh[6]={-35.81, -0.35, 35.11, 71.07, 108.23, 140.49};
01305 
01306   double philow00[6]={-60.27, -32.91, 0.35, 33.61, 67.37, 102.23}; //Ring0 L0
01307   double phihgh00[6]={-33.61, -0.35, 32.91, 66.67, 101.53, 129.49};
01308 
01309   double philow01[6]={-64.67, -34.91, 0.35, 35.61, 71.37, 108.33}; //Ring0 L1
01310   double phihgh01[6]={-35.61, -0.35, 34.91, 70.67, 107.63, 138.19};
01311 
01312 
01313   iring = -10;
01314 
01315   double tmpdy =  std::abs(yhor1);
01316   for (int i=0; i<netabin; i++) {
01317     if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
01318       ietaho = i+1; 
01319       float tmp1 = fabs(tmpdy-etalow[i]);
01320       float tmp2 = fabs(tmpdy-etahgh[i]);
01321  
01322       localyhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01323 
01324       if (i<4) iring =0;
01325       if (i>=4 && i<10) iring=1;
01326       if (i>=10 && i<netabin) iring=2;
01327       break;
01328     }
01329   }
01330 
01331   int tmpphi = 0;
01332   int tmpphi0 = 0;
01333 
01334   if (ietaho >4) { //Ring 1 and 2
01335     for (int i=0; i<6; i++) {
01336       if (xhor1 >philow[i] && xhor1 <phihgh[i]) { 
01337         tmpphi=i+1; 
01338         float tmp1 = fabs(xhor1-philow[i]);
01339         float tmp2 = fabs(xhor1-phihgh[i]);
01340         localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01341         break;
01342       }
01343     }
01344   } else {  //Ring 0
01345     for (int i=0; i<6; i++) {
01346       if (xhor1 >philow01[i] && xhor1 <phihgh01[i]) { 
01347         tmpphi=i+1; 
01348         float tmp1 = fabs(xhor1-philow01[i]);
01349         float tmp2 = fabs(xhor1-phihgh01[i]);
01350         localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01351         break;
01352       }
01353     }
01354 
01355     for (int i=0; i<6; i++) {
01356       if (xhor0 >philow00[i] && xhor0 <phihgh00[i]) { 
01357         tmpphi0=i+1; 
01358         float tmp1 = fabs(xhor0-philow00[i]);
01359         float tmp2 = fabs(xhor0-phihgh00[i]);
01360         localxhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01361         if (tmpphi !=tmpphi0) localxhor0 +=10000.;
01362         break;
01363       }
01364     }
01365 
01366     double tmpdy =  std::abs(yhor0);
01367     for (int i=0; i<4; i++) {
01368       if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
01369         float tmp1 = fabs(tmpdy-etalow[i]);
01370         float tmp2 = fabs(tmpdy-etahgh[i]);
01371         localyhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01372         if (i+1 != ietaho)  localyhor0 +=10000.;
01373         break;
01374       }
01375     }
01376   }
01377 
01378   if (tmpphi!=0) { 
01379     iphiho = 6*iphisect -2 + tmpphi;
01380     if (iphiho <=0) iphiho +=nphimx;
01381     if (iphiho >nphimx) iphiho -=nphimx;
01382   }
01383 
01384   //  isect2 = 15*iring+iphisect+1;
01385 
01386   if (yhor1 <0) { 
01387     if (std::abs(ietaho) >netabin) { //Initialised with 50
01388       ietaho +=1; 
01389     } else {
01390       ietaho *=-1; 
01391     }
01392     //    isect2 *=-1; 
01393     iring *=-1;
01394   } 
01395 }
01396 
01397 FreeTrajectoryState AlCaHOCalibProducer::getFreeTrajectoryState( const reco::Track& tk, const MagneticField* field, int iiner, bool dir)
01398 {
01399 
01400   if (iiner ==0) {
01401     GlobalPoint gpos( tk.outerX(), tk.outerY(), tk.outerZ());
01402     GlobalVector gmom( tk.outerPx(), tk.outerPy(), tk.outerPz());
01403     if (dir) gmom *=-1.;
01404     GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
01405     CurvilinearTrajectoryError err( tk.extra()->outerStateCovariance());
01406     return FreeTrajectoryState( par, err);
01407   } else {
01408     GlobalPoint gpos( tk.innerPosition().X(), tk.innerPosition().Y(), tk.innerPosition().Z());
01409     GlobalVector gmom( tk.innerMomentum().X(), tk.innerMomentum().Y(), tk.innerMomentum().Z());
01410     if (dir) gmom *=-1.;
01411     GlobalTrajectoryParameters par( gpos, -gmom, tk.charge(), field);
01412     CurvilinearTrajectoryError err( tk.extra()->innerStateCovariance());
01413     return FreeTrajectoryState( par, err);
01414   }
01415 
01416 }
01417 
01418 #include "FWCore/Framework/interface/MakerMacros.h"
01419 
01420 //define this as a plug-in
01421 DEFINE_FWK_MODULE(AlCaHOCalibProducer);
01422