CMS 3D CMS Logo

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.18 2009/02/09 16:31:28 kodolova 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/ParameterSet/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 "PhysicsTools/UtilAlgos/interface/TFileService.h"
00129 
00130 #include "DataFormats/L1GlobalTrigger/interface/L1GlobalTriggerReadoutRecord.h"
00131 #include "DataFormats/Common/interface/TriggerResults.h"
00132 #include "FWCore/Framework/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(const edm::EventSetup&) ;
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       m_shape = (*conditions_).getHcalShape();
00427 
00428       for (int i=0; i<netamx; i++) {
00429         for (int j=0; j<nphimx; j++) {
00430           for (int k=0; k<ncidmx; k++) {
00431             pedestal[i][j][k]=0.0;
00432           }
00433         }
00434       }     
00435     }
00436   }
00437 
00438   //  if (m_hotime && m_digiInput) {
00439   if (m_digiInput) {
00440     if (irunold !=irun) {
00441       nRuns++;
00442       for (int i =-netabin+1; i<=netabin-1; i++) {
00443         if (i==0) continue;
00444         int tmpeta1 =  (i>0) ? i -1 : -i +14; 
00445         if (tmpeta1 <0 || tmpeta1 >netamx) continue;
00446         for (int j=0; j<nphimx; j++) {
00447           
00448           HcalDetId id(HcalOuter, i, j+1, 4);
00449           calibped = conditions_->getHcalCalibrations(id);
00450           
00451           for (int k =0; k<ncidmx-1; k++) {
00452             pedestal[tmpeta1][j][k] = calibped.pedestal(k); // pedm->getValue(k);
00453             pedestal[tmpeta1][j][ncidmx-1] += (1./(ncidmx-1))*pedestal[tmpeta1][j][k];
00454           }
00455           
00456           if (m_hotime) {
00457             for (int k =0; k<ncidmx; k++) {
00458               libhoped->Fill(nphimx*ncidmx*tmpeta1 + ncidmx*j + k, pedestal[tmpeta1][j][k]);
00459             }
00460             for (int k =0; k<nchnmx; k++) {
00461               libhoped1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*j + k, pedestal[tmpeta1][j][min(k,ncidmx-1)]);
00462             }
00463           }
00464 
00465         }
00466       }
00467     }
00468   }
00469 
00470   //  Nevents++;
00471   irunold = irun;
00472 
00473   //GMA  if (Nevents%500==1) 
00474   //GMA  cout <<"AlCaHOCalibProducer Processing event # "<<Nevents<<" "<<Npass<<" "<<Noccu<<" "<<irun<<" "<<iEvent.id().event()<<endl;
00475 
00476   std::auto_ptr<HOCalibVariableCollection> hostore (new HOCalibVariableCollection);
00477 
00478   edm::Handle<HODigiCollection> ho;   
00479   bool isHO = true;
00480   
00481   edm::Handle<HBHEDigiCollection> hbhe; 
00482   bool isHBHE = true;
00483 
00484   if (m_digiInput) {
00485     try {
00486       iEvent.getByLabel(hoLabel_,ho);
00487     } catch ( cms::Exception &iEvent ) { isHO = false; } 
00488     
00489     try {
00490       iEvent.getByLabel(hbheLabel_,hbhe);
00491     } catch ( cms::Exception &iEvent ) { isHBHE = false; }  
00492   }
00493   
00494   if (m_hotime && m_digiInput) {
00495     if (isHO && (*ho).size()>0) {
00496       for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
00497         HcalDetId id =(*j).id();
00498         int tmpeta= id.ieta();
00499         int tmpphi= id.iphi();
00500         m_coder = (*conditions_).getHcalCoder(id);
00501         float tmpdata[nchnmx];
00502         int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
00503         for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00504           tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00505           allhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
00506           Nallhotime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.);
00507         }
00508       }
00509     }
00510     if (isHBHE && (*hbhe).size()>0) {
00511       for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
00512         HcalDetId id =(*j).id();
00513         int tmpeta= id.ieta();
00514         int tmpphi= id.iphi();
00515         int tmpdepth =id.depth();
00516         m_coder = (*conditions_).getHcalCoder(id);
00517         int tmpeta1 =  (tmpeta>0) ? tmpeta -15 : -tmpeta + 1; 
00518         if (tmpdepth==1) tmpeta1 =  (tmpeta>0) ? tmpeta -1 : -tmpeta +29;  
00519         for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00520           float signal = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00521           if (tmpdepth==1) { 
00522             allhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00523             Nallhb1->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);
00524             hb1pedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);}
00525           if (tmpdepth==2) { 
00526             allhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00527             Nallhb2->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
00528           if (tmpdepth==3) { 
00529             allhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, signal);
00530             Nallhb3->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1);}
00531         }
00532       }
00533     }
00534   }
00535 
00536   double pival = acos(-1.);
00537   
00538   Handle<reco::TrackCollection> cosmicmuon;
00539   bool isMuon = true;
00540   try {
00541     iEvent.getByLabel(muonTags_, cosmicmuon);
00542     
00543   } catch ( cms::Exception &iEvent ) { isMuon = false; } 
00544   
00545   if (isMuon && cosmicmuon->size()>0) { 
00546     
00547     int l1trg = 0;
00548     int hlttr = 0;
00549     
00550     int ntrgpas_gm[ntrgp_gm]={0,0,0,0,0,0,0,0,0,0};
00551     
00552     //L1 trigger
00553     Handle<L1GlobalTriggerReadoutRecord> L1GTRR;
00554     bool isL1Trig=true;
00555     try{   
00556       iEvent.getByLabel(l1Label_,L1GTRR);  //gtDigis
00557     } catch (cms::Exception &iEvent) { isL1Trig=false;}
00558     
00559     if (isL1Trig && L1GTRR.isValid()) {
00560       const unsigned int n(L1GTRR->decisionWord().size());
00561       const bool accept(L1GTRR->decision());
00562       if (accept) {
00563         for (unsigned int i=0; i!=n && i<32; ++i) {
00564           //    for (unsigned int i=0; i!=n ; ++i) {
00565           int il1trg = (L1GTRR->decisionWord()[i]) ? 1 : 0;
00566           if (il1trg>0 && i<32) l1trg +=int(pow(2., double(i%32))*il1trg);
00567         }
00568       }
00569     }// else { return;}
00570     
00571     //HLT 
00572 
00573     Handle<edm::TriggerResults> trigRes;
00574     bool isTrig=true;
00575 
00576     try{
00577       iEvent.getByLabel(hltLabel_, trigRes);
00578     } catch (cms::Exception &iEvent) { isTrig=false;}
00579     if (isTrig) {
00580       unsigned int size = trigRes->size();
00581       edm::TriggerNames triggerNames(*trigRes);
00582       
00583       // loop over all paths, get trigger decision
00584       for(unsigned i = 0; i != size && i<32; ++i) {
00585         std::string name = triggerNames.triggerName(i);
00586         fired[name] = trigRes->accept(i);
00587         int ihlt =  trigRes->accept(i);
00588         if (m_hotime){ 
00589           if (ihlt >0 && i < (int)ntrgp_gm) { ntrgpas_gm[i] = 1;}
00590         }
00591         if (i<32 && ihlt>0) hlttr += int(pow(2., double(i%32))*ihlt);
00592       }
00593     }
00594 
00595     int Noccu_old = Noccu;
00596     
00597     for(reco::TrackCollection::const_iterator ncosm = cosmicmuon->begin();
00598         ncosm != cosmicmuon->end();  ++ncosm) {
00599       
00600       if ((*ncosm).ndof() < 15) continue;
00601       if ((*ncosm).normalizedChi2() >30.0) continue;
00602 
00603       HOCalibVariables tmpHOCalib;
00604       
00605       tmpHOCalib.trig1 = l1trg;
00606       tmpHOCalib.trig2 = hlttr;    
00607       
00608       int charge = ncosm->charge();  
00609       
00610       double innerr = (*ncosm).innerPosition().Perp2();
00611       double outerr = (*ncosm).outerPosition().Perp2();
00612       int iiner = (innerr <outerr) ? 1 : 0;
00613       
00614       //---------------------------------------------------
00615       //             in_to_out  Dir         in_to_out  Dir
00616       //   StandAlone ^         ^     Cosmic    ^    |
00617       //              |         |               |    v
00618       //---------------------------------------------------Y=0
00619       //   StandAlone |         |     Cosmic    ^    |
00620       //              v         v               |    v
00621       //----------------------------------------------------
00622       
00623       double posx, posy, posz;
00624       double momx, momy, momz;
00625       
00626       if (iiner==1) {
00627         posx = (*ncosm).innerPosition().X();
00628         posy = (*ncosm).innerPosition().Y();
00629         posz = (*ncosm).innerPosition().Z();
00630         
00631         momx = (*ncosm).innerMomentum().X();
00632         momy = (*ncosm).innerMomentum().Y();
00633         momz = (*ncosm).innerMomentum().Z();
00634         
00635       } else {
00636         posx = (*ncosm).outerPosition().X();
00637         posy = (*ncosm).outerPosition().Y();
00638         posz = (*ncosm).outerPosition().Z();
00639         
00640         momx = (*ncosm).outerMomentum().X();
00641         momy = (*ncosm).outerMomentum().Y();
00642         momz = (*ncosm).outerMomentum().Z();
00643       }
00644       
00645       
00646       PositionType trkpos(posx, posy, posz);
00647       
00648       Hep3Vector tmpmuon3v(posx, posy, posz);
00649       Hep3Vector tmpmuondir(momx, momy, momz);
00650       
00651       bool samedir = (tmpmuon3v.dot(tmpmuondir) >0) ? true : false;
00652       for (int i=0; i<3; i++) {tmpHOCalib.caloen[i] = 0.0;}
00653       int inearbymuon = 0;
00654       for(reco::TrackCollection::const_iterator ncosmcor = cosmicmuon->begin();
00655           ncosmcor != cosmicmuon->end();  ++ncosmcor) {
00656         if (ncosmcor==ncosm) continue;
00657         
00658         Hep3Vector tmpmuon3vcor;
00659         Hep3Vector tmpmom3v;
00660         if (iiner==1) {
00661           tmpmuon3vcor = Hep3Vector((*ncosmcor).innerPosition().X(),(*ncosmcor).innerPosition().Y(),(*ncosmcor).innerPosition().Z());
00662           tmpmom3v = Hep3Vector((*ncosmcor).innerMomentum().X(),(*ncosmcor).innerMomentum().Y(),(*ncosmcor).innerMomentum().Z());
00663         } else {
00664           tmpmuon3vcor = Hep3Vector((*ncosmcor).outerPosition().X(),(*ncosmcor).outerPosition().Y(),(*ncosmcor).outerPosition().Z());
00665           tmpmom3v = Hep3Vector((*ncosmcor).outerMomentum().X(),(*ncosmcor).outerMomentum().Y(),(*ncosmcor).outerMomentum().Z());       
00666           
00667         }
00668         if (tmpmom3v.mag()<0.2 || (*ncosmcor).ndof()<5) continue;
00669         
00670         double angle = tmpmuon3v.angle(tmpmuon3vcor);
00671         if (angle < 7.5*pival/180.) {inearbymuon=1;} //  break;}
00672 
00673         if (muonTags_.label() =="cosmicMuons") {
00674           if (angle <7.5*pival/180.) { tmpHOCalib.caloen[0] +=1.;}
00675           if (angle <15.0*pival/180.) { tmpHOCalib.caloen[1] +=1.;}
00676           if (angle <35.0*pival/180.) { tmpHOCalib.caloen[2] +=1.;}
00677         }
00678       }
00679       
00680       localxhor0 = localyhor0 = 20000;  //GM for 22OCT07 data
00681       
00682       if (muonTags_.label() =="standAloneMuons") {
00683         
00684         Handle<CaloTowerCollection> calotower;
00685         bool isCaloTower = true;
00686         
00687         
00688         try {
00689           iEvent.getByLabel(towerLabel_, calotower);
00690         } catch ( cms::Exception &iEvent ) { isCaloTower = false; } 
00691         if (isCaloTower) {
00692           for (CaloTowerCollection::const_iterator calt = calotower->begin();
00693                calt !=calotower->end(); calt++) {
00694             //CMSSW_2_1_x       const math::XYZVector towermom = (*calt).momentum();
00695             double ith = (*calt).momentum().theta();
00696             double iph = (*calt).momentum().phi();
00697 
00698             Hep3Vector calo3v(sin(ith)*cos(iph), sin(ith)*sin(iph), cos(ith));
00699             
00700             double angle = tmpmuon3v.angle(calo3v);
00701             
00702             if (angle < 7.5*pival/180.) {tmpHOCalib.caloen[0] += calt->emEnergy()+calt->hadEnergy();}
00703             if (angle < 15*pival/180.) {tmpHOCalib.caloen[1] += calt->emEnergy()+calt->hadEnergy();}
00704             if (angle < 35*pival/180.) {tmpHOCalib.caloen[2] += calt->emEnergy()+calt->hadEnergy();}
00705           }
00706         } 
00707         
00708       }
00709       if (tmpHOCalib.caloen[0] >10.0) continue;
00710       
00711       GlobalPoint glbpt(posx, posy, posz);
00712       
00713       double mom = sqrt(momx*momx + momy*momy +momz*momz);
00714       
00715       momx /= mom;
00716       momy /= mom;
00717       momz /= mom;
00718       
00719       DirectionType trkdir(momx, momy, momz);
00720       
00721       tmpHOCalib.trkdr = (*ncosm).d0();
00722       tmpHOCalib.trkdz = (*ncosm).dz();
00723       
00724       tmpHOCalib.nmuon = cosmicmuon->size();
00725       tmpHOCalib.trkvx = glbpt.x();
00726       tmpHOCalib.trkvy = glbpt.y();
00727       tmpHOCalib.trkvz = glbpt.z();
00728       tmpHOCalib.trkmm = mom*charge;
00729       tmpHOCalib.trkth = trkdir.theta();
00730       tmpHOCalib.trkph = trkdir.phi();
00731       
00732       tmpHOCalib.ndof  = (inearbymuon ==0) ? (int)(*ncosm).ndof() : -(int)(*ncosm).ndof();
00733       tmpHOCalib.chisq = (*ncosm).normalizedChi2(); // max(1.,tmpHOCalib.ndof);
00734       tmpHOCalib.therr = 0.;
00735       tmpHOCalib.pherr = 0.;
00736       
00737       if (iiner==1) {
00738         reco::TrackBase::CovarianceMatrix innercov = (*ncosm).innerStateCovariance();
00739         tmpHOCalib.therr = innercov(1,1); //thetaError();
00740         tmpHOCalib.pherr = innercov(2,2); //phi0Error();
00741       } else {
00742         reco::TrackBase::CovarianceMatrix outercov = (*ncosm).outerStateCovariance();
00743         tmpHOCalib.therr = outercov(1,1); //thetaError();
00744         tmpHOCalib.pherr = outercov(2,2); //phi0Error();
00745       }
00746       
00747       ESHandle<MagneticField> theMagField;
00748       iSetup.get<IdealMagneticFieldRecord>().get(theMagField );     
00749       GlobalVector magfld = theMagField->inInverseGeV(glbpt);
00750 
00751 
00752       SteppingHelixPropagator myHelix(&*theMagField,anyDirection);
00753       myHelix.setMaterialMode(false);
00754       myHelix.applyRadX0Correction(true);
00755 
00756       double phiho = trkpos.phi();
00757       if (phiho<0) phiho +=2*pival;
00758       
00759       int iphisect_dt=int(6*(phiho+pival/18.)/pival); //for u 18/12/06
00760       if (iphisect_dt>=12) iphisect_dt=0;
00761 
00762       int iphisect = -1;
00763 
00764       int ipath = 0;
00765       for (int kl = 0; kl<=2; kl++) {
00766 
00767         int iphisecttmp = (kl<2) ? iphisect_dt + kl : iphisect_dt - 1;
00768         if (iphisecttmp <0) iphisecttmp = 11;
00769         if (iphisecttmp >=12) iphisecttmp = 0;
00770         
00771         double phipos = iphisecttmp*pival/6.;
00772         double phirot = phipos;
00773         
00774         GlobalVector xLocal(-sin(phirot), cos(phirot), 0.);
00775         
00776         GlobalVector yLocal(0., 0., 1.);
00777         GlobalVector zLocal = xLocal.cross(yLocal).unit();
00778         //    GlobalVector zLocal(cos(phirot), sin(phirot), 0.0); 
00779         
00780 
00781         FreeTrajectoryState freetrajectorystate_ = getFreeTrajectoryState(*ncosm,&(*theMagField), iiner, samedir);
00782         
00783         Surface::RotationType rot(xLocal, yLocal, zLocal);
00784         
00785         for (int ik=1; ik>=0; ik--) { //propagate track in two HO layers
00786           
00787           double radial = 407.0;
00788           if (ik==0) radial = 382.0;
00789 
00790           Surface::PositionType pos(radial*cos(phipos), radial*sin(phipos), 0.);
00791           PlaneBuilder::ReturnType aPlane = PlaneBuilder().plane(pos,rot);
00792 
00793           Surface* aPlane2 = new Plane(pos,rot);
00794 
00795           SteppingHelixStateInfo steppingHelixstateinfo_ = myHelix.propagate(freetrajectorystate_, (*aPlane2));
00796 
00797           if (steppingHelixstateinfo_.isValid()) {
00798 
00799             GlobalVector hotrkpos2(steppingHelixstateinfo_.position().x(), steppingHelixstateinfo_.position().y(), steppingHelixstateinfo_.position().z());
00800             Hep3Vector hotrkdir2(steppingHelixstateinfo_.momentum().x(), steppingHelixstateinfo_.momentum().y(),steppingHelixstateinfo_.momentum().z());
00801             
00802             LocalVector lclvt0 = (*aPlane).toLocal(hotrkpos2);
00803             
00804             double xx = lclvt0.x();
00805             double yy = lclvt0.y();
00806             
00807             if (ik ==1) {
00808               if ((abs(yy) < 130 && xx >-64.7 && xx <138.2)
00809                   ||(abs(yy) > 130 && abs(yy) <700 && xx >-76.3 && xx <140.5)) {
00810                 ipath = 1;  //Only look for tracks which as hits in layer 1
00811                 iphisect = iphisecttmp;
00812               }
00813             }
00814             
00815             if (iphisect != iphisecttmp) continue; //Look for ring-0 only when ring1 is accepted for that sector
00816             
00817             switch (ik) 
00818               {
00819               case 0 : 
00820                 xhor0 = xx; //lclvt0.x();
00821                 yhor0 = yy; //lclvt0.y();
00822                 break;
00823               case 1 :
00824                 xhor1 = xx; //lclvt0.x();
00825                 yhor1 = yy; //lclvt0.y();
00826                 
00827                 tmpHOCalib.hoang = Hep3Vector(zLocal.x(),zLocal.y(),zLocal.z()).dot(hotrkdir2.unit());
00828                 break;
00829               default : break;
00830               }
00831           } else {
00832             break;
00833           }
00834         }
00835         if (ipath) break;
00836       }
00837       if (ipath) { //If muon crossed HO laeyrs
00838         
00839         int ietaho = 50;
00840         int iphiho = -1;
00841         
00842         for (int i=0; i<9; i++) {tmpHOCalib.hosig[i]=-100.0;}
00843         for (int i=0; i<18; i++) {tmpHOCalib.hocorsig[i]=-100.0;}
00844         for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
00845         tmpHOCalib.hocro = -100;
00846         tmpHOCalib.htime = -1000;
00847         
00848         int isect = 0;
00849 
00850         findHOEtaPhi(iphisect, ietaho, iphiho);
00851         
00852         if (ietaho !=0 && iphiho !=0 && abs(iring)<=2) { //Muon passed through a tower
00853           isect = 100*abs(ietaho+30)+abs(iphiho);
00854           if (abs(ietaho) >=netabin || iphiho<0) isect *=-1; //Not extrapolated to any tower
00855           if (abs(ietaho) >=netabin) isect -=1000000;  //not matched with eta
00856           if (iphiho<0)        isect -=2000000; //not matched with phi
00857           tmpHOCalib.isect = isect;
00858           
00859           tmpHOCalib.hodx = localxhor1;
00860           tmpHOCalib.hody = localyhor1;      
00861           
00862           if (iring==0) {
00863             tmpHOCalib.hocorsig[8] = localxhor0;
00864             tmpHOCalib.hocorsig[9] = localyhor0;
00865           }
00866           
00867           int etamn=-4;
00868           int etamx=4;
00869           if (iring==1) {etamn=5; etamx = 10;}
00870           if (iring==2) {etamn=11; etamx = 16;}
00871           if (iring==-1){etamn=-10; etamx = -5;}
00872           if (iring==-2){etamn=-16; etamx = -11;}
00873           
00874           int phimn = 1;
00875           int phimx = 2;
00876           if (iring ==0) {
00877             phimx =2*int((iphiho+1)/2.);
00878             phimn = phimx - 1;
00879           } else {
00880             phimn = 3*int((iphiho+1)/3.) - 1; 
00881             phimx = phimn + 2;
00882           }
00883           
00884           if (phimn <1) phimn += nphimx;
00885           if (phimx >72) phimx -= nphimx;
00886           
00887           int sigstr = m_startTS; // 5;
00888           int sigend = m_endTS; // 8;
00889           
00890           //      if (iphiho <=nphimx/2) { //GMA310508
00891           //        sigstr -=1; //GMA300608
00892           //        sigend -=1;
00893           //      }
00894           
00895           if (m_hbinfo) {
00896             for (int i=0; i<9; i++) {tmpHOCalib.hbhesig[i]=-100.0;}
00897             
00898             if (m_digiInput) {
00899               if (isHBHE && (*hbhe).size() >0) {
00900                 for (HBHEDigiCollection::const_iterator j=(*hbhe).begin(); j!=(*hbhe).end(); j++){
00901                   //              const HBHEDataFrame digi = (const HBHEDataFrame)(*j);
00902                   //              HcalDetId id =digi.id();
00903                   HcalDetId id =(*j).id();
00904                   int tmpeta= id.ieta();
00905                   int tmpphi= id.iphi();
00906                   m_coder = (*conditions_).getHcalCoder(id);
00907                   calibped = conditions_->getHcalCalibrations(id);
00908                   
00909                   int deta = tmpeta-ietaho;
00910                   if (tmpeta==-1 && ietaho== 1) deta = -1;
00911                   if (tmpeta== 1 && ietaho==-1) deta =  1;
00912                   int dphi = tmpphi-iphiho;
00913                   if (phimn >phimx) {
00914                     if (dphi==71) dphi=-1;
00915                     if (dphi==-71) dphi=1;
00916                   }
00917                   
00918                   int ipass2 = (abs(deta) <=1 && abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
00919                   
00920                   if (ipass2 ==0 ) continue;
00921                   
00922                   float tmpdata[nchnmx];
00923                   for (int i=0; i<(*j).size() && i<nchnmx; i++) {
00924                     tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
00925                   }
00926                   
00927                   float signal = 0;
00928                   for (int i=1; i<(*j).size() && i<=8; i++) {
00929                     signal += tmpdata[i] - calibped.pedestal((*j).sample(i).capid());; 
00930                   }
00931                   
00932                   if (ipass2 == 1) {
00933                     if (3*(deta+1)+dphi+1<9)  tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
00934                   }
00935                 }
00936               }
00937               
00938             } else {
00939               
00940               edm::Handle<HBHERecHitCollection> hbheht;// iEvent.getByType(hbheht);
00941               bool isHBHEht = true;
00942               try {
00943                 iEvent.getByLabel(hbheLabel_,hbheht);
00944               } catch ( cms::Exception &iEvent ) { isHBHEht = false; }  
00945               
00946               if (isHBHEht && (*hbheht).size()>0) {
00947                 if(!(*hbheht).size()) throw (int)(*hbheht).size();
00948                 
00949                 for (HBHERecHitCollection::const_iterator j=(*hbheht).begin(); j!=(*hbheht).end(); j++){
00950                   //              const HBHERecHit hbhehtrec = (const HBHERecHit)(*j);
00951                   //              HcalDetId id =hbhehtrec.id();
00952                   HcalDetId id =(*j).id();
00953                   int tmpeta= id.ieta();
00954                   int tmpphi= id.iphi();
00955                   
00956                   int deta = tmpeta-ietaho;
00957                   if (tmpeta==-1 && ietaho== 1) deta = -1;
00958                   if (tmpeta== 1 && ietaho==-1) deta =  1;
00959                   int dphi = tmpphi-iphiho;
00960                   if (phimn >phimx) {
00961                     if (dphi==71) dphi=-1;
00962                     if (dphi==-71) dphi=1;
00963                   }
00964                   
00965                   int ipass2 = (abs(deta) <=1 && abs(dphi)<=1) ? 1 : 0; //NEED correction in full CMS detector
00966                   if ( ipass2 ==0 ) continue;
00967                   
00968                   float signal = (*j).energy();
00969                   
00970                   if (3*(deta+1)+dphi+1<9)  tmpHOCalib.hbhesig[3*(deta+1)+dphi+1] = signal;
00971                 }
00972               }
00973               
00974             } //else m_digilevel
00975             
00976           } //m_hbinfo #endif
00977           
00978           if (m_digiInput) {
00979             if (isHO && (*ho).size()>0) {
00980               int isFilled[netamx*nphimx]; 
00981               for (int j=0; j<netamx*nphimx; j++) {isFilled[j]=0;}
00982               
00983               double sumEt = 0;
00984               double sumE  = 0;
00985               
00986               for (HODigiCollection::const_iterator j=(*ho).begin(); j!=(*ho).end(); j++){
00987                 //              const HODataFrame digi = (const HODataFrame)(*j);
00988                 //              HcalDetId id =digi.id();
00989 
00990                 HcalDetId id =(*j).id();                
00991                 int tmpeta= id.ieta();
00992                 int tmpphi= id.iphi();
00993                 m_coder = (*conditions_).getHcalCoder(id);
00994                 
00995                 int ipass1 =0;
00996                 if (tmpeta >=etamn && tmpeta <=etamx) {
00997                   if (phimn < phimx) {
00998                     ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
00999                   } else {
01000                     ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
01001                   }
01002                 }
01003                 
01004                 int deta = tmpeta-ietaho;
01005                 if (tmpeta==-1 && ietaho== 1) deta = -1;
01006                 if (tmpeta== 1 && ietaho==-1) deta =  1;
01007                 
01008                 int dphi = tmpphi -iphiho;
01009                 if (phimn>phimx) {
01010                   if (dphi==71) dphi=-1;
01011                   if (dphi==-71) dphi=1;
01012                 }
01013                 
01014                 int ipass2 = (abs(deta) <=1 && abs(dphi)<=1) ? 1 : 0;
01015                 
01016                 int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
01017                 
01018                 float tmpdata[nchnmx]={0,0,0,0,0,0,0,0,0,0};
01019                 float sigvall[nsigpk]={0,0,0,0,0,0,0};
01020              
01021                 for (int i=0; i<(*j).size() && i<nchnmx; i++) {
01022                   tmpdata[i] = m_coder->charge(*m_shape,(*j).sample(i).adc(),(*j).sample(i).capid());
01023                   if (deta==0 && dphi==0) { 
01024                     double tmpE = tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01025                     if (tmpE >0) {
01026                       sumEt += i*tmpE;
01027                       sumE  += tmpE;
01028                     }
01029                     if (m_hotime) {
01030                       //calculate signals in 4 time slices, 0-3,.. 6-9
01031                       if (i>=7-nsigpk) {
01032                         for (int ncap=0; ncap<nsigpk; ncap++) {
01033                           if (i-ncap >= nstrbn && i-ncap <= nstrbn+3) { 
01034                             sigvall[ncap] +=tmpdata[i];
01035                           }
01036                         }
01037                         }
01038                       if (i==(*j).size()-1) {
01039                         float mxled=-1;
01040                         int imxled = 0;
01041                         for (int ij=0; ij<nsigpk; ij++) {
01042                           if (sigvall[ij] > mxled) {mxled = sigvall[ij]; imxled=ij;}
01043                         }
01044                         double pedx = 0.0;
01045                         for (int ij=0; ij<4; ij++) {
01046                           pedx +=pedestal[tmpeta1][tmpphi-1][ij];
01047                         }
01048                         if (mxled-pedx >2 && mxled-pedx <20 ) {
01049                           hopeak[ntrgp_gm]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
01050                           for (int jk=0; jk<ntrgp_gm; jk++) {
01051                             if (ntrgpas_gm[jk]>0) {
01052                               hopeak[jk]->Fill(nphimx*tmpeta1 + tmpphi-1, imxled+nstrbn);
01053                             }
01054                           }
01055                           if (tmpdata[5]+tmpdata[6] >1) {
01056                             horatio->Fill(nphimx*tmpeta1 + tmpphi-1, (tmpdata[5]-tmpdata[6])/(tmpdata[5]+tmpdata[6]));
01057                           }
01058                           for (int ij=0; ij<(*j).size() && ij<nchnmx; ij++) {
01059                             hotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
01060                             Nhotime[ntrgp_gm]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
01061                             for (int jk=0; jk<ntrgp_gm; jk++) {
01062                               if (ntrgpas_gm[jk]>0) {
01063                                 hotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, tmpdata[ij]);
01064                                 Nhotime[jk]->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + ij, 1.);
01065                               }
01066                             }
01067                           }
01068                         }
01069                       }
01070                     }
01071                   }
01072                 }
01073 
01074                 if (abs(tmpeta) <=15 && deta==0 && dphi ==0) { 
01075                   float signal = 0;
01076                   int icnt = 0;
01077                   for (int i =0; i<nchnmx && i< (*j).size(); i++) {
01078                     if (i >=sigstr && i<=sigend) continue;
01079                     signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01080                     if (++icnt >=4) break;
01081                   }
01082                   tmpHOCalib.hocro = signal;
01083                 }
01084                 
01085                 if (m_hotime) { 
01086                   if (ipass1 ==0 && ipass2 ==0 && cosmicmuon->size()<=2) {
01087                     if (abs(ietaho) <=netabin && iphiho >0) {
01088                       if ((iphiho >=1 && iphiho<=nphimx/2 && tmpphi >=1 && tmpphi <=nphimx/2) ||
01089                           (iphiho >nphimx/2 && iphiho<=nphimx && tmpphi >nphimx/2 && tmpphi <=nphimx)) {
01090                         if (isFilled[nphimx*tmpeta1+tmpphi-1]==0) {
01091                           isFilled[nphimx*tmpeta1+tmpphi-1]=1;
01092                           for (int i=0; i<(*j).size() && i<nchnmx; i++) {
01093                             hopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
01094                             Nhopedtime->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, 1.); 
01095                             hopedpr->Fill(nphimx*nchnmx*tmpeta1 + nchnmx*(tmpphi-1) + i, tmpdata[i]);
01096                           }
01097                         } //isFilled
01098                       }
01099                     }
01100                   }
01101                 }
01102                 
01103                 if (ipass1 ==0 && ipass2 ==0 ) continue;
01104                 
01105                 float signal = 0;
01106                 for (int i=sigstr; i<(*j).size() && i<=sigend; i++) {
01107                   signal += tmpdata[i] - pedestal[tmpeta1][tmpphi-1][(*j).sample(i).capid()];
01108                 }
01109                 if (signal <-100 || signal >100000) signal = -100;
01110                 if (m_hotime) {
01111                   if (signal >-100 && Noccu == Noccu_old) {
01112                     for (int i=0; i<5; i++) {
01113                       if (signal >(i+2)*m_sigma) {
01114                         ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
01115                       }
01116                     }
01117                   }
01118                 }
01119 
01120                 if (ipass1 ==0 && ipass2 ==0 ) continue;
01121                 
01122                 if (ipass1 ==1) {
01123                   int tmpdph = tmpphi-phimn;
01124                   if (tmpdph<0) tmpdph = 2;  //only case of iphi==1, where phimn=71
01125                   
01126                   int ilog = 2*(tmpeta-etamn)+tmpdph;
01127                   if (iring !=0) { 
01128                     if (iring >0) {
01129                       ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
01130                     } else {
01131                       ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
01132                     }
01133                   }
01134                   if (ilog>-1 && ilog<18) { 
01135                     tmpHOCalib.hocorsig[ilog] = signal;
01136                   }
01137                 }             
01138                 
01139                 if (ipass2 ==1) {
01140                   if (3*(deta+1)+dphi+1<9) tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
01141                 }
01142                 
01143                 /*
01144                 // Possibility to store pedestal by shifting phi tower by 6
01145                 // But, due to missing tower at +-5, we do not have always proper
01146                 // statistics and also in pedestal subtracted data, we do not have
01147                 // signal in that tower
01148                 // 
01149                 if (deta==0 && dphi ==0) {
01150                   int crphi = tmpphi + 6;
01151                   if (crphi >72) crphi -=72;
01152                   
01153                   for (HODigiCollection::const_iterator jcr=(*ho).begin(); jcr!=(*ho).end(); jcr++){
01154                   //                const HODataFrame (*jcr) = (const HODataFrame)(*jcr);
01155                   //                HcalDetId idcr =(*jcr).id();
01156                   HcalDetId id =(*jcr).id();
01157                     int etacr= idcr.ieta();
01158                     int phicr= idcr.iphi();
01159                     m_coder = (*conditions_).getHcalCoder(idcr);
01160                     
01161                     if (tmpeta==etacr && crphi ==phicr) {
01162                       
01163                       float tmpdatacr[nchnmx];
01164                       for (int i=0; i<(*jcr).size() && i<nchnmx; i++) {
01165                         tmpdatacr[i] = m_coder->charge(*m_shape,(*jcr).sample(i).adc(),(*jcr).sample(i).capid());
01166                       }
01167                     }
01168                     }
01169                     }
01170                 */
01171                 
01172             }
01173             tmpHOCalib.htime = sumEt/max(sumE,1.e-6);
01174           } 
01175         } else {
01176           edm::Handle<HORecHitCollection> hoht;
01177           bool isHOht = true;
01178           try {
01179             iEvent.getByLabel(hoLabel_,hoht);
01180           } catch ( cms::Exception &iEvent ) { isHOht = false; } 
01181           
01182           if (isHOht && (*hoht).size()>0) {
01183             for (HORecHitCollection::const_iterator j=(*hoht).begin(); j!=(*hoht).end(); j++){
01184               //                const HORecHit hohtrec = (const HORecHit)(*j);
01185               //                HcalDetId id =hohtrec.id();
01186               HcalDetId id =(*j).id();
01187               int tmpeta= id.ieta();
01188               int tmpphi= id.iphi();
01189 
01190               int ipass1 =0;
01191               if (tmpeta >=etamn && tmpeta <=etamx) {
01192                 if (phimn < phimx) {
01193                   ipass1 = (tmpphi >=phimn && tmpphi <=phimx ) ? 1 : 0;
01194                 } else {
01195                   ipass1 = (tmpphi==71 || tmpphi ==72 || tmpphi==1) ? 1 : 0;
01196                 }
01197               }
01198               
01199               int deta = tmpeta-ietaho;
01200               if (tmpeta==-1 && ietaho== 1) deta = -1;
01201               if (tmpeta== 1 && ietaho==-1) deta =  1;
01202               
01203               int dphi = tmpphi -iphiho;
01204               if (phimn>phimx) {
01205                 if (dphi==71) dphi=-1;
01206                 if (dphi==-71) dphi=1;
01207               }
01208               
01209               float signal = (*j).energy();
01210               if (m_hotime) {
01211                 int tmpeta1 = (tmpeta>0) ? tmpeta -1 : -tmpeta +14; 
01212                 if (signal >-100 && Noccu == Noccu_old) {
01213                   for (int i=0; i<5; i++) {
01214                     if (signal >(i+2)*m_sigma) {
01215                       ho_occupency[i]->Fill(nphimx*tmpeta1+tmpphi-1);
01216                     }
01217                   }
01218                 }
01219               }
01220               
01221               int ipass2 = (abs(deta) <=1 && abs(dphi)<=1) ? 1 : 0;
01222               
01223               if (ipass1 ==0 && ipass2 ==0 ) continue;
01224               
01225               if (ipass1 ==1) {
01226                 int tmpdph = tmpphi-phimn;
01227                 if (tmpdph<0) tmpdph = 2;  //only case of iphi==1, where phimn=71
01228                   
01229                 int ilog = 2*(tmpeta-etamn)+tmpdph;
01230                 if (iring !=0) { 
01231                   if (iring >0) {
01232                     ilog = 3*(tmpeta-etamn)+tmpdph; //Again CMS correction
01233                   } else {
01234                     ilog = 3*(etamx-tmpeta)+tmpdph; //Again CMS correction
01235                   }
01236                 }
01237                 if (ilog>-1 && ilog<18) {
01238                   tmpHOCalib.hocorsig[ilog] = signal;
01239                 }
01240               }       
01241               
01242               if (ipass2 ==1) {
01243                 
01244                 if (3*(deta+1)+dphi+1<9) {
01245                   tmpHOCalib.hosig[3*(deta+1)+dphi+1] = signal; //Again CMS azimuthal near phi 1&72
01246                 }
01247               }
01248               
01249               if (deta==0 && dphi ==0) {
01250                 tmpHOCalib.htime = (*j).time();
01251                 int crphi = tmpphi + 6;
01252                 if (crphi >72) crphi -=72;
01253                 
01254                 for (HORecHitCollection::const_iterator jcr=(*hoht).begin(); jcr!=(*hoht).end(); jcr++){
01255                   const HORecHit reccr = (const HORecHit)(*jcr);
01256                   HcalDetId idcr =reccr.id();
01257                   int etacr= idcr.ieta();
01258                   int phicr= idcr.iphi();
01259                   if (tmpeta==etacr && crphi ==phicr) {
01260                     
01261                     tmpHOCalib.hocro = reccr.energy();
01262                     
01263                   }
01264                 }
01265               }
01266             }
01267           } 
01268         }
01269         
01270         //GMA     Npass++;
01271         if (Noccu == Noccu_old) Noccu++;
01272         hostore->push_back(tmpHOCalib); 
01273         
01274       }
01275     }
01276     
01277     } 
01278   } 
01279 
01280   if (hostore->size()>0) iEvent.put(hostore, "HOCalibVariableCollection");
01281   
01282 }
01283 
01284 // ------------ method called once each job just before starting event loop  ------------
01285 void 
01286 AlCaHOCalibProducer::beginJob(const edm::EventSetup& iSetup)
01287 {
01288   //GMA  Nevents = 0;
01289   //GMA  Npass = 0;
01290   //GMA  Noccu = 0;
01291 
01292   irunold = -1;
01293   nRuns = 0;
01294   //  edm::ESHandle<MagneticField> bField;
01295   //  iSetup.get<IdealMagneticFieldRecord>().get(bField);
01296   //  stepProp  = new SteppingHelixPropagator(&*bField,anyDirection);
01297   //  stepProp->setMaterialMode(false);
01298   //  stepProp->applyRadX0Correction(true);
01299   
01300   for (int i=0; i<netamx; i++) {
01301     for (int j=0; j<nphimx; j++) {
01302       for (int k=0; k<ncidmx; k++) {
01303         pedestal[i][j][k]=0.0;
01304       }
01305     }
01306   }
01307 
01308 
01309 }
01310 
01311 // ------------ method called once each job just after ending the event loop  ------------
01312 void 
01313 AlCaHOCalibProducer::endJob() {
01314 
01315 
01316 }
01317 
01318 void AlCaHOCalibProducer::findHOEtaPhi(int iphisect, int& ietaho, int& iphiho) {
01319   
01320   //18/12/06 : use only position, not angle phi
01321 
01322 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};
01323 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};
01324 
01325   double philow[6]={-76.27, -35.11, 0.35, 35.81, 71.77, 108.93};  //Ring+/-1 & 2
01326   double phihgh[6]={-35.81, -0.35, 35.11, 71.07, 108.23, 140.49};
01327 
01328   double philow00[6]={-60.27, -32.91, 0.35, 33.61, 67.37, 102.23}; //Ring0 L0
01329   double phihgh00[6]={-33.61, -0.35, 32.91, 66.67, 101.53, 129.49};
01330 
01331   double philow01[6]={-64.67, -34.91, 0.35, 35.61, 71.37, 108.33}; //Ring0 L1
01332   double phihgh01[6]={-35.61, -0.35, 34.91, 70.67, 107.63, 138.19};
01333 
01334 
01335   iring = -10;
01336 
01337   double tmpdy =  abs(yhor1);
01338   for (int i=0; i<netabin; i++) {
01339     if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
01340       ietaho = i+1; 
01341       float tmp1 = fabs(tmpdy-etalow[i]);
01342       float tmp2 = fabs(tmpdy-etahgh[i]);
01343  
01344       localyhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01345 
01346       if (i<4) iring =0;
01347       if (i>=4 && i<10) iring=1;
01348       if (i>=10 && i<netabin) iring=2;
01349       break;
01350     }
01351   }
01352 
01353   int tmpphi = 0;
01354   int tmpphi0 = 0;
01355 
01356   if (ietaho >4) { //Ring 1 and 2
01357     for (int i=0; i<6; i++) {
01358       if (xhor1 >philow[i] && xhor1 <phihgh[i]) { 
01359         tmpphi=i+1; 
01360         float tmp1 = fabs(xhor1-philow[i]);
01361         float tmp2 = fabs(xhor1-phihgh[i]);
01362         localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01363         break;
01364       }
01365     }
01366   } else {  //Ring 0
01367     for (int i=0; i<6; i++) {
01368       if (xhor1 >philow01[i] && xhor1 <phihgh01[i]) { 
01369         tmpphi=i+1; 
01370         float tmp1 = fabs(xhor1-philow01[i]);
01371         float tmp2 = fabs(xhor1-phihgh01[i]);
01372         localxhor1 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01373         break;
01374       }
01375     }
01376 
01377     for (int i=0; i<6; i++) {
01378       if (xhor0 >philow00[i] && xhor0 <phihgh00[i]) { 
01379         tmpphi0=i+1; 
01380         float tmp1 = fabs(xhor0-philow00[i]);
01381         float tmp2 = fabs(xhor0-phihgh00[i]);
01382         localxhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01383         if (tmpphi !=tmpphi0) localxhor0 +=10000.;
01384         break;
01385       }
01386     }
01387 
01388     double tmpdy =  abs(yhor0);
01389     for (int i=0; i<4; i++) {
01390       if (tmpdy >etalow[i] && tmpdy <etahgh[i]) {
01391         float tmp1 = fabs(tmpdy-etalow[i]);
01392         float tmp2 = fabs(tmpdy-etahgh[i]);
01393         localyhor0 = (tmp1 < tmp2) ? -tmp1 : tmp2;
01394         if (i+1 != ietaho)  localyhor0 +=10000.;
01395         break;
01396       }
01397     }
01398   }
01399 
01400   if (tmpphi!=0) { 
01401     iphiho = 6*iphisect -2 + tmpphi;
01402     if (iphiho <=0) iphiho +=nphimx;
01403     if (iphiho >nphimx) iphiho -=nphimx;
01404   }
01405 
01406   //  isect2 = 15*iring+iphisect+1;
01407 
01408   if (yhor1 <0) { 
01409     if (abs(ietaho) >netabin) { //Initialised with 50
01410       ietaho +=1; 
01411     } else {
01412       ietaho *=-1; 
01413     }
01414     //    isect2 *=-1; 
01415     iring *=-1;
01416   } 
01417 }
01418 
01419 FreeTrajectoryState AlCaHOCalibProducer::getFreeTrajectoryState( const reco::Track& tk, const MagneticField* field, int iiner, bool dir)
01420 {
01421 
01422   if (iiner ==0) {
01423     GlobalPoint gpos( tk.outerX(), tk.outerY(), tk.outerZ());
01424     GlobalVector gmom( tk.outerPx(), tk.outerPy(), tk.outerPz());
01425     if (dir) gmom *=-1.;
01426     GlobalTrajectoryParameters par( gpos, gmom, tk.charge(), field);
01427     CurvilinearTrajectoryError err( tk.extra()->outerStateCovariance());
01428     return FreeTrajectoryState( par, err);
01429   } else {
01430     GlobalPoint gpos( tk.innerPosition().X(), tk.innerPosition().Y(), tk.innerPosition().Z());
01431     GlobalVector gmom( tk.innerMomentum().X(), tk.innerMomentum().Y(), tk.innerMomentum().Z());
01432     if (dir) gmom *=-1.;
01433     GlobalTrajectoryParameters par( gpos, -gmom, tk.charge(), field);
01434     CurvilinearTrajectoryError err( tk.extra()->innerStateCovariance());
01435     return FreeTrajectoryState( par, err);
01436   }
01437 
01438 }
01439 
01440 #include "FWCore/Framework/interface/MakerMacros.h"
01441 
01442 //define this as a plug-in
01443 DEFINE_FWK_MODULE(AlCaHOCalibProducer);
01444 

Generated on Tue Jun 9 17:25:34 2009 for CMSSW by  doxygen 1.5.4