CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
calcTopMass Class Reference
Inheritance diagram for calcTopMass:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &iEvent, const edm::EventSetup &iSetup) override
 
 calcTopMass (const edm::ParameterSet &)
 
 ~calcTopMass () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ESProxyIndex const * esGetTokenIndices (edm::Transition iTrans) const
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
void updateLookup (eventsetup::ESRecordsToProxyIndices const &)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Attributes

float bMass
 
float cMass
 
TH1F * hMassCorFl0
 
TH1F * hMassCorFlM
 
TH1F * hMassCorMix
 
TH1F * hMassNoCorr
 
std::string m_bT_CorrectorName
 
std::string m_cT_CorrectorName
 
std::string m_qT_CorrectorName
 
std::string m_tT_CorrectorName
 
float qMass
 
edm::InputTag sourceByRefer_
 
edm::InputTag sourcePartons_
 
edm::Handle< reco::JetMatchedPartonsCollectiontheJetPartonMatch
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 55 of file CalcTopMassExample.cc.

Constructor & Destructor Documentation

calcTopMass::calcTopMass ( const edm::ParameterSet iConfig)
explicit

Definition at line 86 of file CalcTopMassExample.cc.

References bMass, cMass, edm::ParameterSet::getParameter(), hMassCorFl0, hMassCorFlM, hMassCorMix, hMassNoCorr, m_bT_CorrectorName, m_cT_CorrectorName, m_qT_CorrectorName, m_tT_CorrectorName, TFileService::make(), qMass, sourceByRefer_, and AlCaHLTBitMon_QueryRunRegistry::string.

86  {
87  sourceByRefer_ = iConfig.getParameter<InputTag>("srcByReference");
88  m_qT_CorrectorName = iConfig.getParameter<std::string>("qTopCorrector");
89  m_cT_CorrectorName = iConfig.getParameter<std::string>("cTopCorrector");
90  m_bT_CorrectorName = iConfig.getParameter<std::string>("bTopCorrector");
91  m_tT_CorrectorName = iConfig.getParameter<std::string>("tTopCorrector");
92 
93  bMass = 4.5;
94  cMass = 1.5;
95  qMass = 0.3;
96 
98  hMassNoCorr = fs->make<TH1F>("hMassNoCorr", "", 100, 100, 300);
99  hMassCorFl0 = fs->make<TH1F>("hMassCorFl0", "", 100, 100, 300);
100  hMassCorFlM = fs->make<TH1F>("hMassCorFlM", "", 100, 100, 300);
101  hMassCorMix = fs->make<TH1F>("hMassCorMix", "", 100, 100, 300);
102 }
T getParameter(std::string const &) const
std::string m_cT_CorrectorName
edm::InputTag sourceByRefer_
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
std::string m_qT_CorrectorName
std::string m_tT_CorrectorName
std::string m_bT_CorrectorName
calcTopMass::~calcTopMass ( )
inlineoverride

Definition at line 58 of file CalcTopMassExample.cc.

References analyze(), and iEvent.

58 {};

Member Function Documentation

void calcTopMass::analyze ( const edm::Event iEvent,
const edm::EventSetup iSetup 
)
override

Definition at line 104 of file CalcTopMassExample.cc.

References edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::begin(), bMass, beam_dqm_sourceclient-live_cfg::cerr, cMass, reco::JetCorrector::correction(), gather_cfg::cout, DEFINE_FWK_MODULE, edm::AssociationVector< KeyRefProd, CVal, KeyRef, SizeType, KeyReferenceHelper >::end(), cppFunctionSkipper::exception, edm::Ref< C, T, F >::get(), edm::Event::getByLabel(), JetCorrector::getJetCorrector(), hMassCorFl0, hMassCorFlM, hMassCorMix, hMassNoCorr, edm::EventBase::id(), edm::Ref< C, T, F >::isNonnull(), dqmiolumiharvest::j, m_bT_CorrectorName, m_cT_CorrectorName, m_qT_CorrectorName, m_tT_CorrectorName, reco::MatchedPartons::physicsDefinitionParton(), qMass, sourceByRefer_, and theJetPartonMatch.

Referenced by ~calcTopMass().

104  {
105  cout << "[calcTopMass] analysing event " << iEvent.id() << endl;
106  try {
108  } catch (std::exception &ce) {
109  cerr << "[calcTopMass] caught std::exception " << ce.what() << endl;
110  return;
111  }
112 
113  // get all correctors from top events
114  const JetCorrector *qTopCorrector = JetCorrector::getJetCorrector(m_qT_CorrectorName, iSetup);
115  const JetCorrector *cTopCorrector = JetCorrector::getJetCorrector(m_cT_CorrectorName, iSetup);
116  const JetCorrector *bTopCorrector = JetCorrector::getJetCorrector(m_bT_CorrectorName, iSetup);
117  const JetCorrector *tTopCorrector = JetCorrector::getJetCorrector(m_tT_CorrectorName, iSetup);
118 
119  TLorentzVector jetsNoCorr[6];
120  TLorentzVector jetsCorFl0[6];
121  TLorentzVector jetsCorFlM[6];
122  TLorentzVector jetsCorMix[6];
123 
124  bool isQuarkFound[6] = {false};
125 
126  for (JetMatchedPartonsCollection::const_iterator j = theJetPartonMatch->begin(); j != theJetPartonMatch->end(); j++) {
127  const math::XYZTLorentzVector theJet = (*j).first.get()->p4();
128  const MatchedPartons aMatch = (*j).second;
129  const GenParticleRef &thePhyDef = aMatch.physicsDefinitionParton();
130 
131  if (thePhyDef.isNonnull()) {
132  int particPDG = thePhyDef.get()->pdgId();
133 
134  if (particPDG == 5) { //b from t
135  double bTcorr = bTopCorrector->correction(theJet);
136  double tTcorr = tTopCorrector->correction(theJet);
137  jetsNoCorr[0].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
138  jetsCorFl0[0].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), 0);
139  jetsCorFlM[0].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), bMass);
140  jetsCorMix[0].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
141  isQuarkFound[0] = true;
142  } else if (particPDG == -5) { //bbar from tbar
143  double bTcorr = bTopCorrector->correction(theJet);
144  double tTcorr = tTopCorrector->correction(theJet);
145  jetsNoCorr[3].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
146  jetsCorFl0[3].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), 0);
147  jetsCorFlM[3].SetPtEtaPhiM(theJet.pt() * bTcorr, theJet.eta(), theJet.phi(), bMass);
148  jetsCorMix[3].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
149  isQuarkFound[3] = true;
150  } else if (particPDG == 2) { //W+ from t
151  double qTcorr = qTopCorrector->correction(theJet);
152  double tTcorr = tTopCorrector->correction(theJet);
153  jetsNoCorr[1].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
154  jetsCorFl0[1].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
155  jetsCorFlM[1].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
156  jetsCorMix[1].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
157  isQuarkFound[1] = true;
158  } else if (particPDG == 4) { //W+ from t
159  double cTcorr = cTopCorrector->correction(theJet);
160  double tTcorr = tTopCorrector->correction(theJet);
161  jetsNoCorr[1].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
162  jetsCorFl0[1].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), 0);
163  jetsCorFlM[1].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), cMass);
164  jetsCorMix[1].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
165  isQuarkFound[1] = true;
166  } else if (particPDG == -1) { //W+ from t
167  double qTcorr = qTopCorrector->correction(theJet);
168  double tTcorr = tTopCorrector->correction(theJet);
169  jetsNoCorr[2].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
170  jetsCorFl0[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
171  jetsCorFlM[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
172  jetsCorMix[2].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
173  isQuarkFound[2] = true;
174  } else if (particPDG == -3) { //W+ from t
175  double qTcorr = qTopCorrector->correction(theJet);
176  double tTcorr = tTopCorrector->correction(theJet);
177  jetsNoCorr[2].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
178  jetsCorFl0[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
179  jetsCorFlM[2].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
180  jetsCorMix[2].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
181  isQuarkFound[2] = true;
182  } else if (particPDG == -2) { //W- from tbar
183  double qTcorr = qTopCorrector->correction(theJet);
184  double tTcorr = tTopCorrector->correction(theJet);
185  jetsNoCorr[4].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
186  jetsCorFl0[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
187  jetsCorFlM[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
188  jetsCorMix[4].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
189  isQuarkFound[4] = true;
190  } else if (particPDG == -4) { //W- from tbar
191  double qTcorr = qTopCorrector->correction(theJet);
192  double tTcorr = tTopCorrector->correction(theJet);
193  jetsNoCorr[4].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
194  jetsCorFl0[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
195  jetsCorFlM[4].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
196  jetsCorMix[4].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
197  isQuarkFound[4] = true;
198  } else if (particPDG == 1) { //W- from tbar
199  double qTcorr = qTopCorrector->correction(theJet);
200  double tTcorr = tTopCorrector->correction(theJet);
201  jetsNoCorr[5].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
202  jetsCorFl0[5].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), 0);
203  jetsCorFlM[5].SetPtEtaPhiM(theJet.pt() * qTcorr, theJet.eta(), theJet.phi(), qMass);
204  jetsCorMix[5].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
205  isQuarkFound[5] = true;
206  } else if (particPDG == 3) { //W- from tbar
207  double cTcorr = cTopCorrector->correction(theJet);
208  double tTcorr = tTopCorrector->correction(theJet);
209  jetsNoCorr[5].SetPtEtaPhiM(theJet.pt(), theJet.eta(), theJet.phi(), theJet.mass());
210  jetsCorFl0[5].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), 0);
211  jetsCorFlM[5].SetPtEtaPhiM(theJet.pt() * cTcorr, theJet.eta(), theJet.phi(), cMass);
212  jetsCorMix[5].SetPtEtaPhiM(theJet.pt() * tTcorr, theJet.eta(), theJet.phi(), 0);
213  isQuarkFound[5] = true;
214  }
215  }
216  }
217 
218  if (isQuarkFound[0] && isQuarkFound[1] && isQuarkFound[2]) {
219  TLorentzVector *theTopPNoCorr = new TLorentzVector(jetsNoCorr[0] + jetsNoCorr[1] + jetsNoCorr[2]);
220  TLorentzVector *theTopPCorFl0 = new TLorentzVector(jetsCorFl0[0] + jetsCorFl0[1] + jetsCorFl0[2]);
221  TLorentzVector *theTopPCorFlM = new TLorentzVector(jetsCorFlM[0] + jetsCorFlM[1] + jetsCorFlM[2]);
222  TLorentzVector *theTopPCorMix = new TLorentzVector(jetsCorMix[0] + jetsCorMix[1] + jetsCorMix[2]);
223  hMassNoCorr->Fill(theTopPNoCorr->M());
224  hMassCorFl0->Fill(theTopPCorFl0->M());
225  hMassCorFlM->Fill(theTopPCorFlM->M());
226  hMassCorMix->Fill(theTopPCorMix->M());
227  }
228 
229  if (isQuarkFound[3] && isQuarkFound[4] && isQuarkFound[5]) {
230  TLorentzVector *theTopMNoCorr = new TLorentzVector(jetsNoCorr[3] + jetsNoCorr[4] + jetsNoCorr[5]);
231  TLorentzVector *theTopMCorFl0 = new TLorentzVector(jetsCorFl0[3] + jetsCorFl0[4] + jetsCorFl0[5]);
232  TLorentzVector *theTopMCorFlM = new TLorentzVector(jetsCorFlM[3] + jetsCorFlM[4] + jetsCorFlM[5]);
233  TLorentzVector *theTopMCorMix = new TLorentzVector(jetsCorMix[3] + jetsCorMix[4] + jetsCorMix[5]);
234  hMassNoCorr->Fill(theTopMNoCorr->M());
235  hMassCorFl0->Fill(theTopMCorFl0->M());
236  hMassCorFlM->Fill(theTopMCorFlM->M());
237  hMassCorMix->Fill(theTopMCorMix->M());
238  }
239 }
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
edm::Handle< reco::JetMatchedPartonsCollection > theJetPartonMatch
std::string m_cT_CorrectorName
edm::InputTag sourceByRefer_
const_iterator end() const
double correction(const LorentzVector &fJet) const
get correction using Jet information only
Definition: JetCorrector.h:47
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
T const * get() const
Returns C++ pointer to the item.
Definition: Ref.h:232
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:488
std::string m_qT_CorrectorName
const GenParticleRef & physicsDefinitionParton() const
static const JetCorrector * getJetCorrector(const std::string &fName, const edm::EventSetup &fSetup)
retrieve corrector from the event setup. troughs exception if something is missing ...
Definition: JetCorrector.cc:48
std::string m_tT_CorrectorName
edm::EventID id() const
Definition: EventBase.h:59
const_iterator begin() const
std::string m_bT_CorrectorName

Member Data Documentation

float calcTopMass::bMass
private

Definition at line 71 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

float calcTopMass::cMass
private

Definition at line 72 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

TH1F* calcTopMass::hMassCorFl0
private

Definition at line 76 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

TH1F* calcTopMass::hMassCorFlM
private

Definition at line 77 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

TH1F* calcTopMass::hMassCorMix
private

Definition at line 78 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

TH1F* calcTopMass::hMassNoCorr
private

Definition at line 75 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

std::string calcTopMass::m_bT_CorrectorName
private

Definition at line 68 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

std::string calcTopMass::m_cT_CorrectorName
private

Definition at line 67 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

std::string calcTopMass::m_qT_CorrectorName
private

Definition at line 66 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

std::string calcTopMass::m_tT_CorrectorName
private

Definition at line 69 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

float calcTopMass::qMass
private

Definition at line 73 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

edm::InputTag calcTopMass::sourceByRefer_
private

Definition at line 63 of file CalcTopMassExample.cc.

Referenced by analyze(), and calcTopMass().

edm::InputTag calcTopMass::sourcePartons_
private

Definition at line 62 of file CalcTopMassExample.cc.

edm::Handle<reco::JetMatchedPartonsCollection> calcTopMass::theJetPartonMatch
private

Definition at line 64 of file CalcTopMassExample.cc.

Referenced by analyze().