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, 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 35 of file SignCaloSpecificAlgo.h.

Constructor & Destructor Documentation

SignCaloSpecificAlgo::SignCaloSpecificAlgo ( )

Definition at line 26 of file SignCaloSpecificAlgo.cc.

References matrix_.

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

Definition at line 32 of file SignCaloSpecificAlgo.cc.

33 {
34 }

Member Function Documentation

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

Definition at line 152 of file SignCaloSpecificAlgo.cc.

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

Referenced by cms::METProducer::produce_CaloMET().

153 {
154 
155  //retreive calo tower information from candidates
156  //start with the first element of the candidate list
157 
158 
159 
160  // 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
161 
162  std::vector<metsig::SigInputObj> signInputVec = makeVectorOutOfCaloTowers(towers, resolutions, noHF, globalThreshold);
163 
164  // now run the significance algorithm.
165 
166  double sign_calo_met_total=0;
167  double sign_calo_met_phi=0;
168  double sign_calo_met_set=0;
169  metsig::significanceAlgo signifalgo;
170  // check the caloMET, if significance was already run continue with the matrix that is stored..
171  signifalgo.addSignifMatrix(matrix_);
172  signifalgo.addObjects(signInputVec);
173  matrix_=signifalgo.getSignifMatrix();
174  significance_ = signifalgo.significance( sign_calo_met_total, sign_calo_met_phi, sign_calo_met_set);
175  // cleanup everything:
176  signInputVec.clear();
177  // and return
178 }
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 44 of file SignCaloSpecificAlgo.h.

References significance_.

Referenced by cms::METProducer::produce_CaloMET().

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

Definition at line 45 of file SignCaloSpecificAlgo.h.

References matrix_.

Referenced by cms::METProducer::produce_CaloMET().

45 {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 51 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().

52 {
53 
54  edm::View<Candidate>::const_iterator towerCand = towers->begin();
55  std::vector<metsig::SigInputObj> signInputVec;
56  //iterate over all CaloTowers and record information
57  for( ; towerCand != towers->end(); towerCand++ ) {
58  const Candidate *candidate = &(*towerCand);
59  if(candidate){
60  const CaloTower * calotower = dynamic_cast<const CaloTower*> (candidate);
61  if(calotower){
62  double sign_tower_et = calotower->et();
63  if(sign_tower_et<globalThreshold)
64  continue;
65  bool wasused=false;
66  double sign_tower_phi = calotower->phi();
67  double sign_tower_sigma_et = 0;
68  double sign_tower_sigma_phi = 0;
69  std::string sign_tower_type = "";
70 
71  bool hadIsDone = false;
72  bool emIsDone = false;
73  int cell = calotower->constituentsSize();
74 
75  while ( --cell >= 0 && (!hadIsDone || !emIsDone) )
76  {
77  DetId id = calotower->constituent( cell );
78  if( !hadIsDone && id.det() == DetId::Hcal )
79  {
80  HcalSubdetector subdet = HcalDetId(id).subdet();
81  if(subdet == HcalBarrel){
82  sign_tower_type = "hadcalotower";
83  sign_tower_et = calotower->hadEt();
84  sign_tower_sigma_et = resolutions.eval(metsig::caloHB,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
85  sign_tower_sigma_phi = resolutions.eval(metsig::caloHB,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
86  }
87  else if(subdet==HcalOuter){
88  sign_tower_type = "hadcalotower";
89  sign_tower_et = calotower->outerEt();
90  sign_tower_sigma_et = resolutions.eval(metsig::caloHO,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
91  sign_tower_sigma_phi = resolutions.eval(metsig::caloHO,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
92  }
93  else if(subdet==HcalEndcap){
94  sign_tower_type = "hadcalotower";
95  sign_tower_et = calotower->hadEt();
96  sign_tower_sigma_et = resolutions.eval(metsig::caloHE,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
97  sign_tower_sigma_phi = resolutions.eval(metsig::caloHE,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
98  }
99  else if(subdet == HcalForward){
100  sign_tower_type = "hadcalotower";
101  sign_tower_et = calotower->et();
102  sign_tower_sigma_et = resolutions.eval(metsig::caloHF,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
103  sign_tower_sigma_phi = resolutions.eval(metsig::caloHF,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
104  }
105  else{
106  edm::LogWarning("SignCaloSpecificAlgo") << " HCAL tower cell not assigned to an HCAL subdetector!!!" << std::endl;
107  }
108  // and book!
109  metsig::SigInputObj temp(sign_tower_type,sign_tower_et,sign_tower_phi,sign_tower_sigma_et,sign_tower_sigma_phi);
110  if(!noHF || subdet !=HcalForward)
111  signInputVec.push_back(temp);
112 
113  wasused=1;
114  hadIsDone = true;
115  }
116  else if( !emIsDone && id.det() == DetId::Ecal )
117  {
118  EcalSubdetector subdet = EcalSubdetector( id.subdetId() );
119 
120  if(subdet == EcalBarrel){
121  sign_tower_type = "emcalotower";
122  sign_tower_et = calotower->emEt();
123  sign_tower_sigma_et = resolutions.eval(metsig::caloEB,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
124  sign_tower_sigma_phi = resolutions.eval(metsig::caloEB,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
125  }
126  else if(subdet == EcalEndcap ){
127  sign_tower_type = "emcalotower";
128  sign_tower_et = calotower->emEt();
129  sign_tower_sigma_et = resolutions.eval(metsig::caloEE,metsig::ET,sign_tower_et,calotower->phi(),calotower->eta());
130  sign_tower_sigma_phi = resolutions.eval(metsig::caloEE,metsig::PHI,sign_tower_et,calotower->phi(),calotower->eta());
131 
132  }
133  else{
134  edm::LogWarning("SignCaloSpecificAlgo") << " ECAL tower cell not assigned to an ECAL subdetector!!!" << std::endl;
135  }
136  metsig::SigInputObj temp(sign_tower_type,sign_tower_et,sign_tower_phi,sign_tower_sigma_et,sign_tower_sigma_phi);
137  signInputVec.push_back(temp);
138  wasused=1;
139  emIsDone = true;
140  }
141  }
142  if(wasused==0)
143  edm::LogWarning("SignCaloSpecificAlgo") << "found non-assigned cell, " << std::endl;
144  }
145  }
146  }
147  return signInputVec;
148 }
boost::indirect_iterator< typename seq_t::const_iterator > const_iterator
Definition: View.h:81
size_t constituentsSize() const
Definition: CaloTower.h:74
HcalSubdetector subdet() const
get the subdetector
Definition: HcalDetId.h:32
DetId constituent(size_t i) const
Definition: CaloTower.h:75
double hadEt() const
Definition: CaloTower.h:85
double outerEt() const
Definition: CaloTower.h:86
double eval(const resolutionType &type, const resolutionFunc &func, const double &et, const double &phi, const double &eta, const double &p) const
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
HcalSubdetector
Definition: HcalAssistant.h:32
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
Definition: DetId.h:20
double et(double vtxZ) const
Definition: CaloTower.h:101
EcalSubdetector
double emEt() const
Definition: CaloTower.h:84
void SignCaloSpecificAlgo::usePreviousSignif ( const std::vector< double > &  values)

Definition at line 36 of file SignCaloSpecificAlgo.cc.

References matrix_.

37 {
38  if(values.size()!=4)
39  return;
40  matrix_(0,0)=values[0];
41  matrix_(0,1)=values[1];
42  matrix_(1,0)=values[2];
43  matrix_(1,1)=values[3];
44  return;
45 }
void SignCaloSpecificAlgo::usePreviousSignif ( const TMatrixD &  matrix)
inline

Member Data Documentation

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

Definition at line 53 of file SignCaloSpecificAlgo.h.

Referenced by calculateBaseCaloMET(), and getSignificance().