CMS 3D CMS Logo

L1GctValidation Class Reference

Description: produces standard plots of Gct output quantities to enable validation of global event quantities in particular. More...

#include <L1Trigger/L1GlobalCaloTrigger/plugins/L1GctValidation.cc>

Inheritance diagram for L1GctValidation:

edm::EDAnalyzer

List of all members.

Public Member Functions

 L1GctValidation (const edm::ParameterSet &)
 ~L1GctValidation ()

Private Member Functions

virtual void analyze (const edm::Event &, const edm::EventSetup &)
virtual void beginJob (const edm::EventSetup &)
virtual void endJob ()

Private Attributes

edm::InputTag m_energy_tag
edm::InputTag m_missHt_tag
TH1F * theHfRing0CountNegativeEta
TH1F * theHfRing0CountPositiveEta
TH1F * theHfRing0EtSumNegativeEta
TH1F * theHfRing0EtSumPositiveEta
TH1F * theHfRing1CountNegativeEta
TH1F * theHfRing1CountPositiveEta
TH1F * theHfRing1EtSumNegativeEta
TH1F * theHfRing1EtSumPositiveEta
std::vector< TH1F * > theJetCounts
TH1F * theMissEtAngle
TH1F * theMissEtInGeV
TH1F * theMissEtInLsb
TH2F * theMissEtVector
TH2F * theMissEtVsMissHt
TH2F * theMissEtVsMissHtAngle
TH1F * theMissHtAngle
TH1F * theMissHtInGeV
TH1F * theMissHtInLsb
TH2F * theMissHtVector
TH1F * theSumEtInGeV
TH1F * theSumEtInLsb
TH1F * theSumHtInGeV
TH1F * theSumHtInLsb


Detailed Description

Description: produces standard plots of Gct output quantities to enable validation of global event quantities in particular.

Definition at line 38 of file L1GctValidation.h.


Constructor & Destructor Documentation

L1GctValidation::L1GctValidation ( const edm::ParameterSet iConfig  )  [explicit]

Definition at line 17 of file L1GctValidation.cc.

00017                                                                :
00018   m_energy_tag(iConfig.getUntrackedParameter<edm::InputTag>("inputTag", edm::InputTag("gctDigis"))),
00019   m_missHt_tag(iConfig.getUntrackedParameter<edm::InputTag>("missHtTag",edm::InputTag("gctDigis:missingHt")))
00020 {
00021 }

L1GctValidation::~L1GctValidation (  ) 

Definition at line 24 of file L1GctValidation.cc.

00025 {
00026  
00027    // do anything here that needs to be done at desctruction time
00028    // (e.g. close files, deallocate resources etc.)
00029 
00030 }


Member Function Documentation

void L1GctValidation::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
) [private, virtual]

Implements edm::EDAnalyzer.

Definition at line 39 of file L1GctValidation.cc.

References funct::cos(), etHad, etMiss, etTot, edm::EventSetup::get(), edm::Event::getByLabel(), m_energy_tag, m_missHt_tag, L1GctJetCounts::MAX_TOTAL_COUNTS, L1GctEtMiss::phi(), funct::sin(), theHfRing0CountNegativeEta, theHfRing0CountPositiveEta, theHfRing0EtSumNegativeEta, theHfRing0EtSumPositiveEta, theHfRing1CountNegativeEta, theHfRing1CountPositiveEta, theHfRing1EtSumNegativeEta, theHfRing1EtSumPositiveEta, theJetCounts, theMissEtAngle, theMissEtInGeV, theMissEtInLsb, theMissEtVector, theMissEtVsMissHt, theMissEtVsMissHtAngle, theMissHtAngle, theMissHtInGeV, theMissHtInLsb, theMissHtVector, theSumEtInGeV, theSumEtInLsb, theSumHtInGeV, and theSumHtInLsb.

00040 {
00041    using namespace edm;
00042 
00043   // Get the scales from the event setup
00044   ESHandle< L1GctJetEtCalibrationFunction > calibFun ;
00045   iSetup.get< L1GctJetCalibFunRcd >().get( calibFun ) ; // which record?
00046   ESHandle< L1CaloEtScale > etScale ;
00047   iSetup.get< L1JetEtScaleRcd >().get( etScale ) ; // which record?
00048 
00049   double lsbForEt = etScale.product()->linearLsb();
00050   double lsbForHt = calibFun.product()->getHtScaleLSB();
00051 
00052   // Get the Gct energy sums from the event
00053   Handle< L1GctEtTotalCollection > sumEtColl ;
00054   iEvent.getByLabel( m_energy_tag, sumEtColl ) ;
00055   Handle< L1GctEtHadCollection >   sumHtColl ;
00056   iEvent.getByLabel( m_energy_tag, sumHtColl ) ;
00057   Handle< L1GctEtMissCollection >  missEtColl ;
00058   iEvent.getByLabel( m_energy_tag, missEtColl ) ;
00059   Handle< L1GctEtMissCollection >  missHtColl ;
00060   iEvent.getByLabel( m_missHt_tag, missHtColl ) ;
00061 
00062   double etTot = 0.0;
00063   for (L1GctEtTotalCollection::const_iterator jbx=sumEtColl->begin(); jbx!=sumEtColl->end(); jbx++) {
00064     if (jbx->bx()==0) { etTot  = static_cast<double>(jbx->et()); }
00065   } 
00066 
00067   double etHad  = 0.0;
00068   for (L1GctEtHadCollection::const_iterator jbx=sumHtColl->begin(); jbx!=sumHtColl->end(); jbx++) {
00069     if (jbx->bx()==0) { etHad  = static_cast<double>(jbx->et()); }
00070   }
00071 
00072   double etMiss = 0.0;
00073   double etMAng = 0.0;
00074   for (L1GctEtMissCollection::const_iterator jbx=missEtColl->begin(); jbx!=missEtColl->end(); jbx++) {
00075     if (jbx->bx()==0) {
00076       etMiss = static_cast<double>(jbx->et());
00077       int phibin = jbx->phi();
00078       if (phibin>=36) phibin -= 72;
00079       double etMPhi = static_cast<double>(phibin);
00080 
00081       etMAng = (etMPhi+0.5)*M_PI/36.;
00082     }
00083   }
00084 
00085   double htMiss = 0.0;
00086   double htMAng = 0.0;
00087   for (L1GctEtMissCollection::const_iterator jbx=missHtColl->begin(); jbx!=missHtColl->end(); jbx++) {
00088     if (jbx->bx()==0) {
00089       htMiss = static_cast<double>(jbx->et());
00090       int phibin = jbx->phi();
00091       if (phibin>=36) phibin -= 72;
00092       double htMPhi = static_cast<double>(phibin);
00093 
00094       htMAng = (htMPhi+0.5)*M_PI/36.;
00095     }
00096   }
00097 
00098   theSumEtInLsb->Fill(etTot);
00099   theSumHtInLsb->Fill(etHad);
00100   theMissEtInLsb->Fill(etMiss);
00101   theMissHtInLsb->Fill(htMiss);
00102   theSumEtInGeV->Fill(etTot*lsbForEt);
00103   theSumHtInGeV->Fill(etHad*lsbForHt);
00104   theMissEtInGeV->Fill(etMiss*lsbForEt);
00105   theMissEtAngle->Fill(etMAng);
00106   theMissEtVector->Fill(etMiss*lsbForEt*cos(etMAng),etMiss*lsbForEt*sin(etMAng));
00107   theMissHtInGeV->Fill(htMiss*lsbForHt*8.);
00108   theMissHtAngle->Fill(htMAng);
00109   theMissHtVector->Fill(htMiss*lsbForHt*cos(htMAng)*8.,htMiss*lsbForHt*sin(htMAng)*8.);
00110 
00111   theMissEtVsMissHt->Fill(etMiss*lsbForEt, htMiss*lsbForHt*8);
00112   theMissEtVsMissHtAngle->Fill(etMAng, htMAng);
00113 
00114   // Get jet counts from the event
00115   Handle< L1GctJetCountsCollection > jetCountColl ;
00116   iEvent.getByLabel( m_energy_tag, jetCountColl ) ;
00117 
00118   for (L1GctJetCountsCollection::const_iterator jbx=jetCountColl->begin(); jbx!=jetCountColl->end(); jbx++) {
00119     if (jbx->bx()==0) {
00120       for (unsigned jc=0; jc<L1GctJetCounts::MAX_TOTAL_COUNTS; jc++) {
00121         theJetCounts.at(jc)->Fill(jbx->count(jc));
00122       }
00123     }
00124 
00125   }
00126 
00127   // Get minbias trigger quantities from HF
00128   Handle<L1GctHFRingEtSumsCollection> HFEtSumsColl;
00129   Handle<L1GctHFBitCountsCollection>  HFCountsColl;
00130   iEvent.getByLabel( m_energy_tag, HFEtSumsColl ) ;
00131   iEvent.getByLabel( m_energy_tag, HFCountsColl ) ;
00132 
00133   for (L1GctHFRingEtSumsCollection::const_iterator es = HFEtSumsColl->begin(); es != HFEtSumsColl->end(); es++) {
00134     if (es->bx()==0) {
00135       theHfRing0EtSumPositiveEta->Fill(es->etSum(0));  
00136       theHfRing0EtSumNegativeEta->Fill(es->etSum(1));  
00137       theHfRing1EtSumPositiveEta->Fill(es->etSum(2));  
00138       theHfRing1EtSumNegativeEta->Fill(es->etSum(3));
00139     }
00140   }  
00141 
00142   for (L1GctHFBitCountsCollection::const_iterator bc = HFCountsColl->begin(); bc != HFCountsColl->end(); bc++) {
00143     if (bc->bx()==0) {
00144       theHfRing0CountPositiveEta->Fill(bc->bitCount(0));  
00145       theHfRing0CountNegativeEta->Fill(bc->bitCount(1));  
00146       theHfRing1CountPositiveEta->Fill(bc->bitCount(2));  
00147       theHfRing1CountNegativeEta->Fill(bc->bitCount(3));
00148     }
00149   }
00150 
00151 }

void L1GctValidation::beginJob ( const edm::EventSetup  )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 155 of file L1GctValidation.cc.

References header, TFileDirectory::make(), L1GctJetCounts::MAX_TOTAL_COUNTS, ss, theHfRing0CountNegativeEta, theHfRing0CountPositiveEta, theHfRing0EtSumNegativeEta, theHfRing0EtSumPositiveEta, theHfRing1CountNegativeEta, theHfRing1CountPositiveEta, theHfRing1EtSumNegativeEta, theHfRing1EtSumPositiveEta, theJetCounts, theMissEtAngle, theMissEtInGeV, theMissEtInLsb, theMissEtVector, theMissEtVsMissHt, theMissEtVsMissHtAngle, theMissHtAngle, theMissHtInGeV, theMissHtInLsb, theMissHtVector, theSumEtInGeV, theSumEtInLsb, theSumHtInGeV, theSumHtInLsb, and indexGen::title.

00156 {
00157   edm::Service<TFileService> fs;
00158 
00159   TFileDirectory dir0 = fs->mkdir("L1GctEtSums");
00160 
00161   theSumEtInLsb   = dir0.make<TH1F>("SumEtInLsb",   "Total Et (GCT units)",
00162                                     128, 0., 2048.); 
00163   theSumHtInLsb   = dir0.make<TH1F>("SumHtInLsb",   "Total Ht (GCT units)",
00164                                     128, 0., 2048.); 
00165   theMissEtInLsb  = dir0.make<TH1F>("MissEtInLsb",  "Missing Et magnitude (GCT units)",
00166                                     128, 0., 1024.); 
00167   theMissHtInLsb  = dir0.make<TH1F>("MissHtInLsb",  "Missing Ht magnitude (GCT units)",
00168                                     64, 0., 63.); 
00169   theSumEtInGeV   = dir0.make<TH1F>("SumEtInGeV",   "Total Et (in GeV)",
00170                                     100, 0., 1000.); 
00171   theSumHtInGeV   = dir0.make<TH1F>("SumHtInGeV",   "Total Ht (in GeV)",
00172                                     100, 0., 1000.); 
00173   theMissEtInGeV  = dir0.make<TH1F>("MissEtInGeV",  "Missing Et magnitude (in GeV)",
00174                                     100, 0., 500.); 
00175   theMissEtAngle  = dir0.make<TH1F>("MissEtAngle",  "Missing Et angle",
00176                                     72, -M_PI, M_PI);
00177   theMissEtVector = dir0.make<TH2F>("MissEtVector", "Missing Ex vs Missing Ey",
00178                                     100, -100., 100., 100, -100., 100.); 
00179   theMissHtInGeV  = dir0.make<TH1F>("MissHtInGeV",  "Missing Ht magnitude (in GeV)",
00180                                     100, 0., 500.); 
00181   theMissHtAngle  = dir0.make<TH1F>("MissHtAngle",  "Missing Ht angle",
00182                                     72, -M_PI, M_PI);
00183   theMissHtVector = dir0.make<TH2F>("MissHtVector", "Missing Hx vs Missing Hy",
00184                                     100, -100., 100., 100, -100., 100.); 
00185   theMissEtVsMissHt = dir0.make<TH2F>("MissEtVsMissHt", "Missing Et vs Missing Ht",
00186                                       100, 0., 500., 100, 0., 500.);
00187   theMissEtVsMissHtAngle = dir0.make<TH2F>("MissEtVsMissHtAngle", "Angle correlation Missing Et vs Missing Ht",
00188                                            72, -M_PI, M_PI, 72, -M_PI, M_PI);
00189 
00190   TFileDirectory dir1 = fs->mkdir("L1GctHfSumsAndJetCounts");
00191 
00192   for (unsigned jc=0; jc<L1GctJetCounts::MAX_TOTAL_COUNTS; jc++) {
00193     std::stringstream ss;
00194     std::string title;
00195     std::string header;
00196     ss << "JetCount#" << jc;
00197     ss >> title;
00198     ss << "Jet Count number " << jc;
00199     if (jc== 6 || jc== 7) { ss << " (Hf tower count)"; }
00200     if (jc== 8 || jc== 9) { ss << " (Hf Et0 sum MSB)"; }
00201     if (jc==10 || jc==11) { ss << " (Hf Et1 sum MSB)"; }
00202     ss >> header;
00203     theJetCounts.push_back(dir1.make<TH1F>(title.c_str(), header.c_str(), 32, 0., 32.));
00204   }
00205 
00206   // !!!+++
00207   // The following needs updating for new Hf sums code
00208   // !!!+++
00209   theHfRing0EtSumPositiveEta = dir1.make<TH1F>("HfRing0EtSumPositiveEta", "Hf Inner Ring0 Et eta+",
00210                                                60, 0., 30.);
00211   theHfRing0EtSumNegativeEta = dir1.make<TH1F>("HfRing0EtSumNegativeEta", "Hf Inner Ring0 Et eta-",
00212                                                60, 0., 30.);
00213   theHfRing1EtSumPositiveEta = dir1.make<TH1F>("HfRing1EtSumPositiveEta", "Hf Inner Ring1 Et eta+",
00214                                                60, 0., 30.);
00215   theHfRing1EtSumNegativeEta = dir1.make<TH1F>("HfRing1EtSumNegativeEta", "Hf Inner Ring1 Et eta-",
00216                                                60, 0., 30.);
00217   theHfRing0CountPositiveEta = dir1.make<TH1F>("HfRing0CountPositiveEta", "Hf Threshold bits Ring0 eta+",
00218                                                20, 0., 20.);
00219   theHfRing0CountNegativeEta = dir1.make<TH1F>("HfRing0CountNegativeEta", "Hf Threshold bits Ring0 eta-",
00220                                                20, 0., 20.);
00221   theHfRing1CountPositiveEta = dir1.make<TH1F>("HfRing1CountPositiveEta", "Hf Threshold bits Ring1 eta+",
00222                                                20, 0., 20.);
00223   theHfRing1CountNegativeEta = dir1.make<TH1F>("HfRing1CountNegativeEta", "Hf Threshold bits Ring1 eta-",
00224                                                20, 0., 20.);
00225 }

void L1GctValidation::endJob ( void   )  [private, virtual]

Reimplemented from edm::EDAnalyzer.

Definition at line 229 of file L1GctValidation.cc.

00229                         {
00230 }


Member Data Documentation

edm::InputTag L1GctValidation::m_energy_tag [private]

Definition at line 51 of file L1GctValidation.h.

Referenced by analyze().

edm::InputTag L1GctValidation::m_missHt_tag [private]

Definition at line 52 of file L1GctValidation.h.

Referenced by analyze().

TH1F* L1GctValidation::theHfRing0CountNegativeEta [private]

Definition at line 76 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing0CountPositiveEta [private]

Definition at line 75 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing0EtSumNegativeEta [private]

Definition at line 72 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing0EtSumPositiveEta [private]

Definition at line 71 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing1CountNegativeEta [private]

Definition at line 78 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing1CountPositiveEta [private]

Definition at line 77 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing1EtSumNegativeEta [private]

Definition at line 74 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theHfRing1EtSumPositiveEta [private]

Definition at line 73 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

std::vector<TH1F*> L1GctValidation::theJetCounts [private]

Definition at line 70 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissEtAngle [private]

Definition at line 61 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissEtInGeV [private]

Definition at line 60 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissEtInLsb [private]

Definition at line 56 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH2F* L1GctValidation::theMissEtVector [private]

Definition at line 62 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH2F* L1GctValidation::theMissEtVsMissHt [private]

Definition at line 67 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH2F* L1GctValidation::theMissEtVsMissHtAngle [private]

Definition at line 68 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissHtAngle [private]

Definition at line 64 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissHtInGeV [private]

Definition at line 63 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theMissHtInLsb [private]

Definition at line 57 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH2F* L1GctValidation::theMissHtVector [private]

Definition at line 65 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theSumEtInGeV [private]

Definition at line 58 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theSumEtInLsb [private]

Definition at line 54 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theSumHtInGeV [private]

Definition at line 59 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().

TH1F* L1GctValidation::theSumHtInLsb [private]

Definition at line 55 of file L1GctValidation.h.

Referenced by analyze(), and beginJob().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:26:39 2009 for CMSSW by  doxygen 1.5.4