CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
SignCaloSpecificAlgo.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: METAlgorithms
4 // Class: SigInputObj
5 //
6 // Original Author: Kyle Story, Freya Blekman (Cornell University)
7 // Created: Fri Apr 18 11:58:33 CEST 2008
8 //
9 //
10 
11 //____________________________________________________________________________||
18 
19 #include <string>
20 
21 using namespace reco;
22 using namespace std;
23 
24 //____________________________________________________________________________||
26  significance_(0.),
27  matrix_(2,2)
28 {
29  matrix_(0,0)=matrix_(1,0)=matrix_(0,1)=matrix_(1,1)=0.;
30 }
32 {
33 }
34 
35 void SignCaloSpecificAlgo::usePreviousSignif(const std::vector<double> &values)
36 {
37  if(values.size()!=4)
38  return;
39  matrix_(0,0)=values[0];
40  matrix_(0,1)=values[1];
41  matrix_(1,0)=values[2];
42  matrix_(1,1)=values[3];
43  return;
44 }
46 //
47 // Convert a list of calo towers to objects that can be passed to the significance algo:
48 
49 std::vector<metsig::SigInputObj>
50 SignCaloSpecificAlgo::makeVectorOutOfCaloTowers(edm::Handle<edm::View<reco::Candidate> > towers, const::metsig::SignAlgoResolutions& resolutions, bool noHF, double globalThreshold)
51 {
52 
53  edm::View<Candidate>::const_iterator towerCand = towers->begin();
54  std::vector<metsig::SigInputObj> signInputVec;
55  //iterate over all CaloTowers and record information
56  for( ; towerCand != towers->end(); towerCand++ ) {
57  const Candidate *candidate = &(*towerCand);
58  if(candidate){
59  const CaloTower * calotower = dynamic_cast<const CaloTower*> (candidate);
60  if(calotower){
61  double sign_tower_et = calotower->et();
62  if(sign_tower_et<globalThreshold)
63  continue;
64  bool wasused=false;
65  double sign_tower_phi = calotower->phi();
66  double sign_tower_sigma_et = 0;
67  double sign_tower_sigma_phi = 0;
68  std::string sign_tower_type = "";
69 
70  bool hadIsDone = false;
71  bool emIsDone = false;
72  int cell = calotower->constituentsSize();
73 
74  while ( --cell >= 0 && (!hadIsDone || !emIsDone) )
75  {
76  DetId id = calotower->constituent( cell );
77  if( !hadIsDone && id.det() == DetId::Hcal )
78  {
79  HcalSubdetector subdet = HcalDetId(id).subdet();
80  if(subdet == HcalBarrel){
81  sign_tower_type = "hadcalotower";
82  sign_tower_et = calotower->hadEt();
83  sign_tower_sigma_et = resolutions.eval(metsig::caloHB,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
84  sign_tower_sigma_phi = resolutions.eval(metsig::caloHB,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
85  }
86  else if(subdet==HcalOuter){
87  sign_tower_type = "hadcalotower";
88  sign_tower_et = calotower->outerEt();
89  sign_tower_sigma_et = resolutions.eval(metsig::caloHO,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
90  sign_tower_sigma_phi = resolutions.eval(metsig::caloHO,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
91  }
92  else if(subdet==HcalEndcap){
93  sign_tower_type = "hadcalotower";
94  sign_tower_et = calotower->hadEt();
95  sign_tower_sigma_et = resolutions.eval(metsig::caloHE,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
96  sign_tower_sigma_phi = resolutions.eval(metsig::caloHE,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
97  }
98  else if(subdet == HcalForward){
99  sign_tower_type = "hadcalotower";
100  sign_tower_et = calotower->et();
101  sign_tower_sigma_et = resolutions.eval(metsig::caloHF,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
102  sign_tower_sigma_phi = resolutions.eval(metsig::caloHF,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
103  }
104  else{
105  edm::LogWarning("SignCaloSpecificAlgo") << " HCAL tower cell not assigned to an HCAL subdetector!!!" << std::endl;
106  }
107  // and book!
108  metsig::SigInputObj temp(sign_tower_type,sign_tower_et,sign_tower_phi,sign_tower_sigma_et,sign_tower_sigma_phi);
109  if(!noHF || subdet !=HcalForward)
110  signInputVec.push_back(temp);
111 
112  wasused=1;
113  hadIsDone = true;
114  }
115  else if( !emIsDone && id.det() == DetId::Ecal )
116  {
117  EcalSubdetector subdet = EcalSubdetector( id.subdetId() );
118 
119  if(subdet == EcalBarrel){
120  sign_tower_type = "emcalotower";
121  sign_tower_et = calotower->emEt();
122  sign_tower_sigma_et = resolutions.eval(metsig::caloEB,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
123  sign_tower_sigma_phi = resolutions.eval(metsig::caloEB,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
124  }
125  else if(subdet == EcalEndcap ){
126  sign_tower_type = "emcalotower";
127  sign_tower_et = calotower->emEt();
128  sign_tower_sigma_et = resolutions.eval(metsig::caloEE,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
129  sign_tower_sigma_phi = resolutions.eval(metsig::caloEE,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
130 
131  }
132  else{
133  edm::LogWarning("SignCaloSpecificAlgo") << " ECAL tower cell not assigned to an ECAL subdetector!!!" << std::endl;
134  }
135  metsig::SigInputObj temp(sign_tower_type,sign_tower_et,sign_tower_phi,sign_tower_sigma_et,sign_tower_sigma_phi);
136  signInputVec.push_back(temp);
137  wasused=1;
138  emIsDone = true;
139  }
140  }
141  if(wasused==0)
142  edm::LogWarning("SignCaloSpecificAlgo") << "found non-assigned cell, " << std::endl;
143  }
144  }
145  }
146  return signInputVec;
147 }
149 //
150 // Basic MET algorithm. gets towers, does sum. Very similar to standard MET.
151 void SignCaloSpecificAlgo::calculateBaseCaloMET(edm::Handle<edm::View<reco::Candidate> > towers, const CommonMETData& met,const metsig::SignAlgoResolutions& resolutions, bool noHF, double globalThreshold)
152 {
153 
154  //retreive calo tower information from candidates
155  //start with the first element of the candidate list
156 
157 
158 
159  // use this container to calculate the significance. SigInputObj are objects that contain both directional and uncertainty information and are used as input to the significance calculation
160 
161  std::vector<metsig::SigInputObj> signInputVec = makeVectorOutOfCaloTowers(towers, resolutions, noHF, globalThreshold);
162 
163  // now run the significance algorithm.
164 
165  double sign_calo_met_total=0;
166  double sign_calo_met_phi=0;
167  double sign_calo_met_set=0;
168  metsig::significanceAlgo signifalgo;
169  // check the caloMET, if significance was already run continue with the matrix that is stored..
170  signifalgo.addSignifMatrix(matrix_);
171  signifalgo.addObjects(signInputVec);
172  matrix_=signifalgo.getSignifMatrix();
173  significance_ = signifalgo.significance( sign_calo_met_total, sign_calo_met_phi, sign_calo_met_set);
174  // cleanup everything:
175  signInputVec.clear();
176  // and return
177 }
178 
179 //-------------------------------------------------------------------------
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
std::vector< metsig::SigInputObj > makeVectorOutOfCaloTowers(edm::Handle< edm::View< reco::Candidate > > towers, const metsig::SignAlgoResolutions &resolutions, bool noHF, double globalthreshold)
tuple met
____________________________________________________________________________||
Definition: CaloMET_cfi.py:7
size_t constituentsSize() const
Definition: CaloTower.h:72
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
DetId constituent(size_t i) const
Definition: CaloTower.h:73
double hadEt() const
Definition: CaloTower.h:83
const double significance(double &met_r, double &met_phi, double &met_set)
double outerEt() const
Definition: CaloTower.h:84
void usePreviousSignif(const std::vector< double > &values)
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
TMatrixD getSignifMatrix() const
HcalSubdetector
Definition: HcalAssistant.h:31
Structure containing data common to all types of MET.
Definition: CommonMETData.h:12
const void addSignifMatrix(const TMatrixD &input)
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
Definition: DetId.h:18
void calculateBaseCaloMET(edm::Handle< edm::View< reco::Candidate > > towers, const CommonMETData &met, const metsig::SignAlgoResolutions &resolutions, bool noHF, double globalthreshold)
double et(double vtxZ) const
Definition: CaloTower.h:99
const_iterator end() const
EcalSubdetector
double emEt() const
Definition: CaloTower.h:82