CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
SignCaloSpecificAlgo Class Reference

#include <RecoMET/METAlgorithms/interface/SignCaloSpecificAlgo.h>

Public Member Functions

void calculateBaseCaloMET (edm::Handle< edm::View< reco::Candidate > > towers, const CommonMETData &met, const metsig::SignAlgoResolutions &resolutions, bool noHF, double globalthreshold)
 
double getSignificance ()
 
TMatrixD getSignificanceMatrix () const
 
 SignCaloSpecificAlgo ()
 
void usePreviousSignif (const std::vector< double > &values)
 
void usePreviousSignif (const TMatrixD &matrix)
 
 ~SignCaloSpecificAlgo ()
 

Private Member Functions

std::vector< metsig::SigInputObjmakeVectorOutOfCaloTowers (edm::Handle< edm::View< reco::Candidate > > towers, const metsig::SignAlgoResolutions &resolutions, bool noHF, double globalthreshold)
 

Private Attributes

TMatrixD matrix_
 
double significance_
 

Detailed Description

Description: <one line="" class="" summary>="">

Implementation: <Notes on="" implementation>="">

Definition at line 34 of file SignCaloSpecificAlgo.h.

Constructor & Destructor Documentation

SignCaloSpecificAlgo::SignCaloSpecificAlgo ( )

Definition at line 25 of file SignCaloSpecificAlgo.cc.

References matrix_.

25  :
26  significance_(0.),
27  matrix_(2,2)
28 {
29  matrix_(0,0)=matrix_(1,0)=matrix_(0,1)=matrix_(1,1)=0.;
30 }
SignCaloSpecificAlgo::~SignCaloSpecificAlgo ( )

Definition at line 31 of file SignCaloSpecificAlgo.cc.

32 {
33 }

Member Function Documentation

void SignCaloSpecificAlgo::calculateBaseCaloMET ( edm::Handle< edm::View< reco::Candidate > >  towers,
const CommonMETData met,
const metsig::SignAlgoResolutions resolutions,
bool  noHF,
double  globalthreshold 
)

Definition at line 151 of file SignCaloSpecificAlgo.cc.

References metsig::significanceAlgo::addObjects(), metsig::significanceAlgo::addSignifMatrix(), metsig::significanceAlgo::getSignifMatrix(), makeVectorOutOfCaloTowers(), matrix_, metsig::significanceAlgo::significance(), and significance_.

Referenced by cms::CaloMETProducer::produce().

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 }
const void addObjects(const std::vector< metsig::SigInputObj > &EventVec)
std::vector< metsig::SigInputObj > makeVectorOutOfCaloTowers(edm::Handle< edm::View< reco::Candidate > > towers, const metsig::SignAlgoResolutions &resolutions, bool noHF, double globalthreshold)
const double significance(double &met_r, double &met_phi, double &met_set)
TMatrixD getSignifMatrix() const
const void addSignifMatrix(const TMatrixD &input)
double SignCaloSpecificAlgo::getSignificance ( )
inline

Definition at line 43 of file SignCaloSpecificAlgo.h.

References significance_.

Referenced by cms::CaloMETProducer::produce().

43 {return significance_;}
TMatrixD SignCaloSpecificAlgo::getSignificanceMatrix ( void  ) const
inline

Definition at line 44 of file SignCaloSpecificAlgo.h.

References matrix_.

Referenced by cms::CaloMETProducer::produce().

44 {return matrix_;}
std::vector< metsig::SigInputObj > SignCaloSpecificAlgo::makeVectorOutOfCaloTowers ( edm::Handle< edm::View< reco::Candidate > >  towers,
const metsig::SignAlgoResolutions resolutions,
bool  noHF,
double  globalthreshold 
)
private

Definition at line 50 of file SignCaloSpecificAlgo.cc.

References metsig::caloEB, metsig::caloEE, metsig::caloHB, metsig::caloHE, metsig::caloHF, metsig::caloHO, CaloTower::constituent(), CaloTower::constituentsSize(), DetId::Ecal, EcalBarrel, EcalEndcap, CaloTower::emEt(), edm::View< T >::end(), metsig::ET, CaloTower::et(), reco::LeafCandidate::eta(), CaloTower::hadEt(), DetId::Hcal, HcalBarrel, HcalEndcap, HcalForward, HcalOuter, CaloTower::outerEt(), metsig::PHI, reco::LeafCandidate::phi(), AlCaHLTBitMon_QueryRunRegistry::string, HcalDetId::subdet(), and groupFilesInBlocks::temp.

Referenced by calculateBaseCaloMET().

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 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
size_t constituentsSize() const
Definition: CaloTower.h:89
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:30
DetId constituent(size_t i) const
Definition: CaloTower.h:90
double hadEt() const
Definition: CaloTower.h:100
virtual float phi() const
momentum azimuthal angle
double outerEt() const
Definition: CaloTower.h:101
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
virtual float eta() const
momentum pseudorapidity
HcalSubdetector
Definition: HcalAssistant.h:31
Definition: DetId.h:18
const_iterator begin() const
double et(double vtxZ) const
Definition: CaloTower.h:116
const_iterator end() const
EcalSubdetector
double emEt() const
Definition: CaloTower.h:99
void SignCaloSpecificAlgo::usePreviousSignif ( const std::vector< double > &  values)

Definition at line 35 of file SignCaloSpecificAlgo.cc.

References matrix_.

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 }
void SignCaloSpecificAlgo::usePreviousSignif ( const TMatrixD &  matrix)
inline

Member Data Documentation

TMatrixD SignCaloSpecificAlgo::matrix_
private
double SignCaloSpecificAlgo::significance_
private

Definition at line 52 of file SignCaloSpecificAlgo.h.

Referenced by calculateBaseCaloMET(), and getSignificance().