CMS 3D CMS Logo

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 ()
 
reco::METCovMatrix getSignificanceMatrix () const
 
 SignCaloSpecificAlgo ()
 
void usePreviousSignif (const std::vector< double > &values)
 
void usePreviousSignif (const reco::METCovMatrix &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

reco::METCovMatrix 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::SignCaloSpecificAlgo ( )

Definition at line 25 of file SignCaloSpecificAlgo.cc.

References matrix_.

25  : significance_(0.) {
26  matrix_(0, 0) = matrix_(1, 0) = matrix_(0, 1) = matrix_(1, 1) = 0.;
27 }
reco::METCovMatrix matrix_

◆ ~SignCaloSpecificAlgo()

SignCaloSpecificAlgo::~SignCaloSpecificAlgo ( )

Definition at line 28 of file SignCaloSpecificAlgo.cc.

28 {}

Member Function Documentation

◆ calculateBaseCaloMET()

void SignCaloSpecificAlgo::calculateBaseCaloMET ( edm::Handle< edm::View< reco::Candidate > >  towers,
const 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(), HLT_2024v12_cff::globalThreshold, makeVectorOutOfCaloTowers(), matrix_, HLT_2024v12_cff::noHF, electronProducer_cfi::resolutions, metsig::significanceAlgo::significance(), significance_, and HLT_2024v12_cff::towers.

156  {
157  //retreive calo tower information from candidates
158  //start with the first element of the candidate list
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)
reco::METCovMatrix matrix_
const void addSignifMatrix(const reco::METCovMatrix &input)
reco::METCovMatrix getSignifMatrix() const

◆ getSignificance()

double SignCaloSpecificAlgo::getSignificance ( )
inline

Definition at line 41 of file SignCaloSpecificAlgo.h.

References significance_.

41 { return significance_; }

◆ getSignificanceMatrix()

reco::METCovMatrix SignCaloSpecificAlgo::getSignificanceMatrix ( void  ) const
inline

Definition at line 42 of file SignCaloSpecificAlgo.h.

References matrix_.

42 { return matrix_; }
reco::METCovMatrix matrix_

◆ makeVectorOutOfCaloTowers()

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 43 of file SignCaloSpecificAlgo.cc.

References metsig::caloEB, metsig::caloEE, metsig::caloHB, metsig::caloHE, metsig::caloHF, metsig::caloHO, ewkMuLumiMonitorDQM_cfi::calotower, DetId::Ecal, EcalBarrel, EcalEndcap, metsig::ET, HLT_2024v12_cff::globalThreshold, DetId::Hcal, HcalBarrel, HcalEndcap, HcalForward, HcalOuter, HLT_2024v12_cff::noHF, metsig::PHI, electronProducer_cfi::resolutions, AlCaHLTBitMon_QueryRunRegistry::string, HcalDetId::subdet(), groupFilesInBlocks::temp, and HLT_2024v12_cff::towers.

Referenced by calculateBaseCaloMET().

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

◆ usePreviousSignif() [1/2]

void SignCaloSpecificAlgo::usePreviousSignif ( const std::vector< double > &  values)

Definition at line 30 of file SignCaloSpecificAlgo.cc.

References matrix_, and contentValuesCheck::values.

30  {
31  if (values.size() != 4)
32  return;
33  matrix_(0, 0) = values[0];
34  matrix_(0, 1) = values[1];
35  matrix_(1, 0) = values[2];
36  matrix_(1, 1) = values[3];
37  return;
38 }
reco::METCovMatrix matrix_

◆ usePreviousSignif() [2/2]

void SignCaloSpecificAlgo::usePreviousSignif ( const reco::METCovMatrix matrix)
inline

Member Data Documentation

◆ matrix_

reco::METCovMatrix SignCaloSpecificAlgo::matrix_
private

◆ significance_

double SignCaloSpecificAlgo::significance_
private

Definition at line 56 of file SignCaloSpecificAlgo.h.

Referenced by calculateBaseCaloMET(), and getSignificance().