ViewVC Help
View File | Revision Log | Show Annotations | Root Listing
root/cvsroot/UserCode/MitAna/DataTree/interface/CaloMet.h
Revision: 1.4
Committed: Wed Mar 18 15:44:32 2009 UTC (16 years, 1 month ago) by loizides
Content type: text/plain
Branch: MAIN
CVS Tags: Mit_011, Mit_010a, Mit_010, Mit_009c, Mit_009b, Mit_009a, Mit_009, Mit_008
Changes since 1.3: +20 -20 lines
Log Message:
Introduced Double32_t [0,0,14] consistently. Updated class descriptions.

File Contents

# User Rev Content
1 ksung 1.1 //--------------------------------------------------------------------------------------------------
2 loizides 1.4 // $Id: CaloMet.h,v 1.3 2009/03/12 15:56:50 bendavid Exp $
3 ksung 1.1 //
4     // CaloMet
5     //
6 loizides 1.4 // Class to hold CaloMet specific information based on caloremetric information.
7 ksung 1.1 //
8 bendavid 1.3 // Authors: C.Loizides
9 ksung 1.1 //--------------------------------------------------------------------------------------------------
10    
11 bendavid 1.3 #ifndef MITANA_DATATREE_CALOMET_H
12     #define MITANA_DATATREE_CALOMET_H
13    
14 ksung 1.1 #include "MitAna/DataTree/interface/Met.h"
15    
16 bendavid 1.3 namespace mithep
17 ksung 1.1 {
18     class CaloMet : public Met
19     {
20     public:
21 bendavid 1.3 CaloMet() :
22     fCaloMetSig(0), fMaxEtInEmTowers(0), fMaxEtInHadTowers(0), fEtFractionHadronic(0),
23     fEmEtFraction(0), fHadEtInHB(0), fHadEtInHO(0), fHadEtInHE(0), fHadEtInHF(0),
24     fEmEtInEB(0), fEmEtInEE(0), fEmEtInHF(0), fCaloSumEtInpHF(0), fCaloSumEtInmHF(0),
25     fCaloMetInpHF(0), fCaloMetInmHF(0), fCaloMetPhiInpHF(0), fCaloMetPhiInmHF(0) {}
26     CaloMet(Double_t mex, Double_t mey) :
27     Met(mex,mey),
28     fCaloMetSig(0), fMaxEtInEmTowers(0), fMaxEtInHadTowers(0), fEtFractionHadronic(0),
29     fEmEtFraction(0), fHadEtInHB(0), fHadEtInHO(0), fHadEtInHE(0), fHadEtInHF(0),
30     fEmEtInEB(0), fEmEtInEE(0), fEmEtInHF(0), fCaloSumEtInpHF(0), fCaloSumEtInmHF(0),
31     fCaloMetInpHF(0), fCaloMetInmHF(0), fCaloMetPhiInpHF(0), fCaloMetPhiInmHF(0) {}
32    
33     Double_t CaloSumEtInpHF() const { return fCaloSumEtInpHF; }
34     Double_t CaloSumEtInmHF() const { return fCaloSumEtInmHF; }
35     Double_t CaloMetInpHF() const { return fCaloMetInpHF; }
36     Double_t CaloMetInmHF() const { return fCaloMetInmHF; }
37     Double_t CaloMetPhiInpHF() const { return fCaloMetPhiInpHF; }
38     Double_t CaloMetPhiInmHF() const { return fCaloMetPhiInmHF; }
39     Double_t EmEtFraction() const { return fEmEtFraction; }
40     Double_t EmEtInEB() const { return fEmEtInEB; }
41     Double_t EmEtInEE() const { return fEmEtInEE; }
42     Double_t EmEtInHF() const { return fEmEtInHF; }
43     Double_t EtFractionHadronic() const { return fEtFractionHadronic; }
44     Double_t HadEtInHB() const { return fHadEtInHB; }
45     Double_t HadEtInHO() const { return fHadEtInHO; }
46     Double_t HadEtInHE() const { return fHadEtInHE; }
47     Double_t HadEtInHF() const { return fHadEtInHF; }
48     Double_t MaxEtInEmTowers() const { return fMaxEtInEmTowers; }
49     Double_t MaxEtInHadTowers() const { return fMaxEtInHadTowers; }
50     Double_t CaloMetSig() const { return fCaloMetSig; }
51     EObjType ObjType() const { return kCaloMet; }
52     void SetCaloSumEtInpHF(Double_t x) { fCaloSumEtInpHF = x; }
53     void SetCaloSumEtInmHF(Double_t x) { fCaloSumEtInmHF = x; }
54     void SetCaloMetInpHF(Double_t x) { fCaloMetInpHF = x; }
55     void SetCaloMetInmHF(Double_t x) { fCaloMetInmHF = x; }
56     void SetCaloMetPhiInpHF(Double_t x) { fCaloMetPhiInpHF = x; }
57     void SetCaloMetPhiInmHF(Double_t x) { fCaloMetPhiInmHF = x; }
58     void SetEmEtFraction(Double_t x) { fEmEtFraction = x; }
59     void SetEtFractionHadronic(Double_t x) { fEtFractionHadronic = x; }
60     void SetHadEtInHB(Double_t x) { fHadEtInHB = x; }
61     void SetHadEtInHO(Double_t x) { fHadEtInHO = x; }
62     void SetHadEtInHE(Double_t x) { fHadEtInHE = x; }
63     void SetHadEtInHF(Double_t x) { fHadEtInHF = x; }
64     void SetEmEtInEB(Double_t x) { fEmEtInEB = x; }
65     void SetEmEtInEE(Double_t x) { fEmEtInEE = x; }
66     void SetEmEtInHF(Double_t x) { fEmEtInHF = x; }
67     void SetMaxEtInEmTowers(Double_t x) { fMaxEtInEmTowers = x; }
68     void SetMaxEtInHadTowers(Double_t x) { fMaxEtInHadTowers = x; }
69     void SetCaloMetSig(Double_t x) { fCaloMetSig = x; }
70 ksung 1.1
71     protected:
72 bendavid 1.3
73 loizides 1.4 Double32_t fCaloMetSig; //[0,0,14]calo met significance
74     Double32_t fMaxEtInEmTowers; //[0,0,14]max ET deposited in ECAL towers
75     Double32_t fMaxEtInHadTowers; //[0,0,14]max ET deposited in HCAL towers
76     Double32_t fEtFractionHadronic; //[0,0,14]event hadronic scaler ET fraction
77     Double32_t fEmEtFraction; //[0,0,14]event em scalar ET fraction
78     Double32_t fHadEtInHB; //[0,0,14]event hadronic scalar ET in HB
79     Double32_t fHadEtInHO; //[0,0,14]event hadronic scalar ET in HO
80     Double32_t fHadEtInHE; //[0,0,14]event hadronic scalar ET in HE
81     Double32_t fHadEtInHF; //[0,0,14]event hadronic scalar ET in HF
82     Double32_t fEmEtInEB; //[0,0,14]event em scalar ET in EB
83     Double32_t fEmEtInEE; //[0,0,14]event em scalar ET in EE
84     Double32_t fEmEtInHF; //[0,0,14]event em scalar ET from HF
85     Double32_t fCaloSumEtInpHF; //[0,0,14]SumET in HF+
86     Double32_t fCaloSumEtInmHF; //[0,0,14]SumET in HF-
87     Double32_t fCaloMetInpHF; //[0,0,14]MET in HF+
88     Double32_t fCaloMetInmHF; //[0,0,14]MET in HF-
89     Double32_t fCaloMetPhiInpHF; //[0,0,14]MET-phi in HF+
90     Double32_t fCaloMetPhiInmHF; //[0,0,14]MET-phi in HF-
91 bendavid 1.3
92     ClassDef(CaloMet, 1) // Missing transverse energy class
93 ksung 1.1 };
94     }
95     #endif