CMS 3D CMS Logo

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

#include <SignCaloSpecificAlgo.h>

Public Types

typedef math::XYZTLorentzVector LorentzVector
 
typedef math::XYZPoint Point
 
typedef std::vector< const
reco::Candidate * > 
TowerCollection
 

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

Definition at line 35 of file SignCaloSpecificAlgo.h.

Member Typedef Documentation

Definition at line 43 of file SignCaloSpecificAlgo.h.

Definition at line 44 of file SignCaloSpecificAlgo.h.

Definition at line 45 of file SignCaloSpecificAlgo.h.

Constructor & Destructor Documentation

SignCaloSpecificAlgo::SignCaloSpecificAlgo ( )

Definition at line 42 of file SignCaloSpecificAlgo.cc.

References matrix_.

42  :
43  significance_(0.),
44  matrix_(2,2)
45 {
46  matrix_(0,0)=matrix_(1,0)=matrix_(0,1)=matrix_(1,1)=0.;
47 }
SignCaloSpecificAlgo::~SignCaloSpecificAlgo ( )

Definition at line 48 of file SignCaloSpecificAlgo.cc.

49 {
50 }

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 168 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().

169 {
170 
171  //retreive calo tower information from candidates
172  //start with the first element of the candidate list
173 
174 
175 
176  // 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
177 
178  std::vector<metsig::SigInputObj> signInputVec = makeVectorOutOfCaloTowers(towers, resolutions, noHF, globalThreshold);
179 
180  // now run the significance algorithm.
181 
182  double sign_calo_met_total=0;
183  double sign_calo_met_phi=0;
184  double sign_calo_met_set=0;
185  metsig::significanceAlgo signifalgo;
186  // check the caloMET, if significance was already run continue with the matrix that is stored..
187  signifalgo.addSignifMatrix(matrix_);
188  signifalgo.addObjects(signInputVec);
189  matrix_=signifalgo.getSignifMatrix();
190  significance_ = signifalgo.significance( sign_calo_met_total, sign_calo_met_phi, sign_calo_met_set);
191  // cleanup everything:
192  signInputVec.clear();
193  // and return
194 }
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 48 of file SignCaloSpecificAlgo.h.

References significance_.

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

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

Definition at line 49 of file SignCaloSpecificAlgo.h.

References matrix_.

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

49 {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 67 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(), HcalDetId::subdet(), and groupFilesInBlocks::temp.

Referenced by calculateBaseCaloMET().

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

Definition at line 52 of file SignCaloSpecificAlgo.cc.

References matrix_.

53 {
54  if(values.size()!=4)
55  return;
56  matrix_(0,0)=values[0];
57  matrix_(0,1)=values[1];
58  matrix_(1,0)=values[2];
59  matrix_(1,1)=values[3];
60  return;
61 }
void SignCaloSpecificAlgo::usePreviousSignif ( const TMatrixD &  matrix)
inline

Member Data Documentation

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

Definition at line 57 of file SignCaloSpecificAlgo.h.

Referenced by calculateBaseCaloMET(), and getSignificance().