CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
FSQDQM Class Reference

#include <FSQDQM.h>

Inheritance diagram for FSQDQM:
DQMEDAnalyzer edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns > edm::one::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 FSQDQM (const edm::ParameterSet &ps)
 
 ~FSQDQM () override
 
- Public Member Functions inherited from DQMEDAnalyzer
void accumulate (edm::Event const &ev, edm::EventSetup const &es) final
 
void beginLuminosityBlock (edm::LuminosityBlock const &lumi, edm::EventSetup const &setup) override
 
void beginRun (edm::Run const &run, edm::EventSetup const &setup) final
 
virtual void dqmBeginRun (edm::Run const &, edm::EventSetup const &)
 
 DQMEDAnalyzer ()
 
 DQMEDAnalyzer (DQMEDAnalyzer const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer &&)=delete
 
void endLuminosityBlock (edm::LuminosityBlock const &, edm::EventSetup const &) override
 
void endLuminosityBlockProduce (edm::LuminosityBlock &lumi, edm::EventSetup const &setup) final
 
void endRun (edm::Run const &run, edm::EventSetup const &setup) override
 
void endRunProduce (edm::Run &run, edm::EventSetup const &setup) override
 
 ~DQMEDAnalyzer () override=default
 
- Public Member Functions inherited from edm::one::EDProducer< edm::Accumulator, edm::EndLuminosityBlockProducer, edm::EndRunProducer, edm::one::WatchLuminosityBlocks, edm::one::WatchRuns >
 EDProducer ()=default
 
SerialTaskQueueglobalLuminosityBlocksQueue () final
 
SerialTaskQueueglobalRunsQueue () final
 
bool hasAbilityToProduceInLumis () const final
 
bool hasAbilityToProduceInRuns () const final
 
bool wantsGlobalLuminosityBlocks () const final
 
bool wantsGlobalRuns () const final
 
- Public Member Functions inherited from edm::one::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
bool wantsStreamLuminosityBlocks () const
 
bool wantsStreamRuns () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
 ~ProducerBase () noexcept(false) 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
 
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)
 
virtual ~EDConsumerBase () noexcept(false)
 

Protected Member Functions

void analyze (edm::Event const &e, edm::EventSetup const &eSetup) override
 
- 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 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)
 

Private Member Functions

void bookHistograms (DQMStore::IBooker &bei, edm::Run const &, edm::EventSetup const &) override
 
void bookHistos (DQMStore *bei)
 

Private Attributes

unsigned int bxNumber_
 
MonitorElementCastorJetMulti
 
MonitorElementCastorJetphi
 
unsigned int eventNumber_
 
MonitorElementh_leadingtrkpt_ntrk_away
 
MonitorElementh_leadingtrkpt_ntrk_towards
 
MonitorElementh_leadingtrkpt_ntrk_transverse
 
MonitorElementh_leadingtrkpt_ptsum_away
 
MonitorElementh_leadingtrkpt_ptsum_towards
 
MonitorElementh_leadingtrkpt_ptsum_transverse
 
MonitorElementh_ntracks
 
MonitorElementh_ntracks_away
 
MonitorElementh_ntracks_towards
 
MonitorElementh_ntracks_transverse
 
MonitorElementh_ptsum_away
 
MonitorElementh_ptsum_towards
 
MonitorElementh_ptsum_transverse
 
MonitorElementh_trkptsum
 
std::vector< int > hltresults
 
std::string labelBS_
 
std::string labelCastorJet_
 
std::string labelPFJet_
 
std::string labelTrack_
 
unsigned int lumiNumber_
 
MonitorElementNPV
 
MonitorElementPFJeteta
 
MonitorElementPFJetMulti
 
MonitorElementPFJetphi
 
MonitorElementPFJetpt
 
MonitorElementPFJetRapidity
 
MonitorElementPV_chi2
 
MonitorElementPV_d0
 
MonitorElementPV_numTrks
 
MonitorElementPV_sumTrks
 
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
 
unsigned int runNumber_
 
edm::EDGetTokenT< reco::BasicJetCollectiontok_castorjet_
 
edm::EDGetTokenT< reco::PFJetCollectiontok_pfjet_
 
edm::EDGetTokenT< reco::TrackCollectiontok_track_
 
MonitorElementTrack_HP_dxyvtx_over_dxyerror
 
MonitorElementTrack_HP_dzvtx_over_dzerr
 
MonitorElementTrack_HP_Eta
 
MonitorElementTrack_HP_Phi
 
MonitorElementTrack_HP_Pt
 
MonitorElementTrack_HP_ptErr_over_pt
 
edm::InputTag vertex_
 

Additional Inherited Members

- Public Types inherited from edm::one::EDProducerBase
typedef EDProducerBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::one::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Attributes inherited from DQMEDAnalyzer
edm::EDPutTokenT< DQMTokenlumiToken_
 
edm::EDPutTokenT< DQMTokenrunToken_
 

Detailed Description

Definition at line 80 of file FSQDQM.h.

Constructor & Destructor Documentation

FSQDQM::FSQDQM ( const edm::ParameterSet ps)

Definition at line 90 of file FSQDQM.cc.

References edm::ParameterSet::getParameter(), and AlCaHLTBitMon_QueryRunRegistry::string.

91 {
92  edm::LogInfo("FSQDQM") << " Creating FSQDQM "
93  << "\n";
94  pvs_ =consumes<edm::View<reco::Vertex> >(
95  iConfig.getParameter<edm::InputTag>("pvs"));
96  labelPFJet_ = iConfig.getParameter<std::string>("LabelPFJet");
97  labelCastorJet_ = iConfig.getParameter<std::string>("LabelCastorJet");
98  labelTrack_ = iConfig.getParameter<std::string>("LabelTrack");
99  tok_pfjet_ = consumes<reco::PFJetCollection>(labelPFJet_);
100  tok_track_ = consumes<reco::TrackCollection>(labelTrack_);
101  tok_castorjet_ = consumes<reco::BasicJetCollection>(labelCastorJet_);
102 
103 }
edm::EDGetTokenT< reco::BasicJetCollection > tok_castorjet_
Definition: FSQDQM.h:100
std::string labelPFJet_
Definition: FSQDQM.h:96
std::string labelTrack_
Definition: FSQDQM.h:96
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
Definition: FSQDQM.h:97
edm::EDGetTokenT< reco::PFJetCollection > tok_pfjet_
Definition: FSQDQM.h:99
std::string labelCastorJet_
Definition: FSQDQM.h:96
edm::EDGetTokenT< reco::TrackCollection > tok_track_
Definition: FSQDQM.h:98
FSQDQM::~FSQDQM ( )
override

Definition at line 106 of file FSQDQM.cc.

107  {
108 
109  edm::LogInfo("FSQDQM") << " Deleting FSQDQM "
110  << "\n";
111  // do anything here that needs to be done at desctruction time
112  // (e.g. close files, deallocate resources etc.)
113 
114  }

Member Function Documentation

void FSQDQM::analyze ( edm::Event const &  e,
edm::EventSetup const &  eSetup 
)
overrideprotectedvirtual

Reimplemented from DQMEDAnalyzer.

Definition at line 162 of file FSQDQM.cc.

References funct::abs(), edm::EventBase::bunchCrossing(), hiPixelPairStep_cff::deltaPhi, edm::EventID::event(), edm::Event::getByToken(), edm::EventBase::id(), reco::Vertex::isFake(), edm::HandleBase::isValid(), gen::k, edm::EventID::luminosityBlock(), reco::Vertex::normalizedChi2(), reco::Vertex::position(), edm::Handle< T >::product(), reco::TrackBase::qualityByName(), RecoPFJets_cff::recoPFJets, edm::EventID::run(), sistrip::runNumber_, mathSSE::sqrt(), reco::Vertex::tracks_begin(), reco::Vertex::tracks_end(), reco::Vertex::tracksSize(), reco::Vertex::x(), reco::Vertex::xError(), reco::Vertex::y(), reco::Vertex::yError(), reco::Vertex::z(), and reco::Vertex::zError().

163  {
164  using namespace edm;
165  using namespace std;
166  using namespace reco;
167 
168  runNumber_ = iEvent.id().run();
169  eventNumber_ = iEvent.id().event();
170  lumiNumber_ = iEvent.id().luminosityBlock();
171  bxNumber_ = iEvent.bunchCrossing();
172 
174  // if(!iEvent.getByToken(pvs_, privtxs)) return;
175  iEvent.getByToken(pvs_, privtxs);
176  double bestvz=-999.9, bestvx=-999.9, bestvy=-999.9;
177  double bestvzError=-999.9, bestvxError=-999.9, bestvyError=-999.9;
178  if(privtxs.isValid()){
179  const reco::Vertex& pvtx=privtxs->front();
180  NPV->Fill(privtxs->size());
181  if(privtxs->begin() !=privtxs->end() && !(pvtx.isFake()) && pvtx.position().Rho() <= 2. && fabs(pvtx.position().z()) <= 24){
182  bestvz = pvtx.z();
183  bestvx = pvtx.x();
184  bestvy = pvtx.y();
185  bestvzError = pvtx.zError();
186  bestvxError = pvtx.xError();
187  bestvyError = pvtx.yError();
188  PV_chi2->Fill(pvtx.normalizedChi2());
189  PV_d0->Fill(sqrt(pvtx.x() * pvtx.x() + pvtx.y() * pvtx.y()));
190  PV_numTrks->Fill(pvtx.tracksSize());
191  double vertex_sumTrks = 0.0;
192  for(reco::Vertex::trackRef_iterator iTrack= pvtx.tracks_begin(); iTrack != pvtx.tracks_end(); iTrack++)
193  {
194  vertex_sumTrks += (*iTrack)->pt();
195  }
196  PV_sumTrks->Fill(vertex_sumTrks);
197  }
198  }
199 
200  std::vector<Jet> recoPFJets;
201  recoPFJets.clear();
202  int nPFCHSJet=0;
203  edm::Handle<PFJetCollection> pfjetchscoll;
204  // if(!iEvent.getByToken(tok_pfjet_, pfjetchscoll)) return;
205  iEvent.getByToken(tok_pfjet_, pfjetchscoll);
206  if(pfjetchscoll.isValid()){
207  const reco::PFJetCollection *pfchsjets = pfjetchscoll.product();
208  reco::PFJetCollection::const_iterator pfjetchsclus = pfchsjets->begin();
209  for(pfjetchsclus = pfchsjets->begin(); pfjetchsclus!= pfchsjets->end() ; ++pfjetchsclus){
210  PFJetpt->Fill( pfjetchsclus->pt());
211  PFJeteta->Fill( pfjetchsclus->eta());
212  PFJetphi->Fill( pfjetchsclus->phi());
213  PFJetRapidity->Fill( pfjetchsclus->rapidity());
214  nPFCHSJet++;
215  }
216  PFJetMulti->Fill( nPFCHSJet);
217  }
218 
219  std::vector<Jet> recoCastorJets;
220  recoCastorJets.clear();
222  // if(!iEvent.getByToken(tok_castorjet_, castorJets)) return;
223  iEvent.getByToken(tok_castorjet_, castorJets);
224  if(castorJets.isValid()){
225  for (unsigned ijet=0; ijet<castorJets->size();ijet++) {
226  recoCastorJets.push_back((*castorJets)[ijet]);
227  }
228  for (unsigned ijet=0; ijet<recoCastorJets.size(); ijet++) {
229  // cout<<recoCastorJets[ijet].pt()<<endl;
230  CastorJetphi->Fill( recoCastorJets[ijet].phi());
231 
232  CastorJetMulti->Fill(recoCastorJets.size());
233  }
234  }
236  // if(!iEvent.getByToken(tok_track_, itracks)) return;
237  iEvent.getByToken(tok_track_, itracks);
238  if(itracks.isValid()){
240  std::vector<TLorentzVector> T_trackRec_P4;
241 
242  int ntracks = 0;
243  int ntracks_towards = 0;
244  int ntracks_transverse = 0;
245  int ntracks_away = 0;
246  float ptsum = 0;
247  float dphi=0;
248  float ptsum_towards = 0;
249  float ptsum_transverse = 0;
250  float ptsum_away = 0;
251 
252  T_trackRec_P4.clear();
253  for(reco::TrackCollection::const_iterator iT = itracks->begin(); iT != itracks->end(); ++iT){
254  if(iT->quality(hiPurity)){
255  math::XYZPoint bestvtx(bestvx,bestvy,bestvz);
256  double dzvtx = iT->dz(bestvtx);
257  double dxyvtx = iT->dxy(bestvtx);
258  double dzerror = sqrt(iT->dzError()*iT->dzError()+bestvzError*bestvzError);
259  double dxyerror = sqrt(iT->d0Error()*iT->d0Error()+bestvxError*bestvyError);
260  if((iT->ptError())/iT->pt() < 0.05 && dzvtx < 3.0 && dxyvtx < 3.0){
261  TLorentzVector trk;
262  trk.SetPtEtaPhiE(iT->pt(),iT->eta(),iT->phi(),iT->p());
263  T_trackRec_P4.push_back(trk);
264  Track_HP_Eta->Fill(iT->eta());
265  Track_HP_Phi->Fill(iT->phi());
266  Track_HP_Pt->Fill(iT->pt());
267  Track_HP_ptErr_over_pt->Fill((iT->ptError())/iT->pt());
268  Track_HP_dzvtx_over_dzerr->Fill(dzvtx/dzerror);
269  Track_HP_dxyvtx_over_dxyerror->Fill(dxyvtx/dxyerror);
270  }
271  }
272  }
273 
274  float highest_pt_track=-999;
275  int index=-999;
276  for(unsigned int k=0;k<T_trackRec_P4.size();k++){
277  if(T_trackRec_P4.at(k).Pt() > highest_pt_track){
278  highest_pt_track=T_trackRec_P4.at(k).Pt();
279  index = k;
280  }
281 
282  }
283  unsigned int finalid=abs(index);
284  // if(T_trackRec_P4.at(index).Pt()!=0.){
285  if(finalid < T_trackRec_P4.size()){
286  // std::sort(T_trackRec_P4.begin(), T_trackRec_P4.end(), SortByPt());
287  for(unsigned int itrk=0;itrk<T_trackRec_P4.size();itrk++){
288  ++ntracks;
289  ptsum= ptsum + T_trackRec_P4.at(itrk).Pt();
290  dphi = deltaPhi(T_trackRec_P4.at(itrk).Phi(),T_trackRec_P4.at(index).Phi());
291  if(fabs(dphi) < 1.05){
292  ++ntracks_towards;
293  ptsum_towards = ptsum_towards + T_trackRec_P4.at(itrk).Pt();}
294  if(fabs(dphi) > 1.05 && fabs(dphi) < 2.09){
295  ++ntracks_transverse;
296  ptsum_transverse = ptsum_transverse + T_trackRec_P4.at(itrk).Pt();}
297  if(fabs(dphi) > 2.09){
298  ++ntracks_away;
299  ptsum_away = ptsum_away + T_trackRec_P4.at(itrk).Pt();}
300  }
301 
302 
303 
304  h_ntracks->Fill(ntracks);
305  h_trkptsum->Fill(ptsum);
306  h_ptsum_towards->Fill(ptsum_towards);
307  h_ptsum_transverse->Fill(ptsum_transverse);
308  h_ptsum_away->Fill(ptsum_away);
309  h_ntracks_towards->Fill(ntracks_towards);
310  h_ntracks_transverse->Fill(ntracks_transverse);
311  h_ntracks_away->Fill(ntracks_away);
312 
313  if(!T_trackRec_P4.empty()){
314  h_leadingtrkpt_ntrk_towards->Fill(T_trackRec_P4.at(index).Pt(),ntracks_towards/8.37);
315  h_leadingtrkpt_ntrk_transverse->Fill(T_trackRec_P4.at(index).Pt(),ntracks_transverse/8.37);
316  h_leadingtrkpt_ntrk_away->Fill(T_trackRec_P4.at(index).Pt(),ntracks_away/8.37);
317  h_leadingtrkpt_ptsum_towards->Fill(T_trackRec_P4.at(index).Pt(),ptsum_towards/8.37);
318  h_leadingtrkpt_ptsum_transverse->Fill(T_trackRec_P4.at(index).Pt(),ptsum_transverse/8.37);
319  h_leadingtrkpt_ptsum_away->Fill(T_trackRec_P4.at(index).Pt(),ptsum_away/8.37);
320  }
321  }
322  }
323  }//analyze
unsigned int lumiNumber_
Definition: FSQDQM.h:103
edm::EDGetTokenT< reco::BasicJetCollection > tok_castorjet_
Definition: FSQDQM.h:100
MonitorElement * h_leadingtrkpt_ntrk_transverse
Definition: FSQDQM.h:137
trackRef_iterator tracks_end() const
last iterator over tracks
Definition: Vertex.cc:81
MonitorElement * PV_d0
Definition: FSQDQM.h:123
double zError() const
error on z
Definition: Vertex.h:123
TrackQuality
track quality
Definition: TrackBase.h:151
MonitorElement * NPV
Definition: FSQDQM.h:121
MonitorElement * PV_numTrks
Definition: FSQDQM.h:124
double y() const
y coordinate
Definition: Vertex.h:113
MonitorElement * Track_HP_dxyvtx_over_dxyerror
Definition: FSQDQM.h:120
MonitorElement * h_ptsum_transverse
Definition: FSQDQM.h:127
MonitorElement * h_ntracks_away
Definition: FSQDQM.h:131
MonitorElement * h_ptsum_towards
Definition: FSQDQM.h:126
MonitorElement * h_ntracks
Definition: FSQDQM.h:133
MonitorElement * PFJetMulti
Definition: FSQDQM.h:113
MonitorElement * h_ntracks_towards
Definition: FSQDQM.h:129
MonitorElement * Track_HP_Phi
Definition: FSQDQM.h:115
const Point & position() const
position
Definition: Vertex.h:109
unsigned int runNumber_
Definition: FSQDQM.h:103
void Fill(long long x)
MonitorElement * PV_sumTrks
Definition: FSQDQM.h:125
MonitorElement * Track_HP_Eta
Definition: FSQDQM.h:116
int iEvent
Definition: GenABIO.cc:230
MonitorElement * PFJetRapidity
Definition: FSQDQM.h:114
MonitorElement * h_trkptsum
Definition: FSQDQM.h:132
MonitorElement * h_leadingtrkpt_ptsum_away
Definition: FSQDQM.h:138
MonitorElement * PFJetphi
Definition: FSQDQM.h:108
T sqrt(T t)
Definition: SSEVec.h:18
edm::EDGetTokenT< edm::View< reco::Vertex > > pvs_
Definition: FSQDQM.h:97
edm::EDGetTokenT< reco::PFJetCollection > tok_pfjet_
Definition: FSQDQM.h:99
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double z() const
z coordinate
Definition: Vertex.h:115
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * h_leadingtrkpt_ptsum_towards
Definition: FSQDQM.h:139
unsigned int eventNumber_
Definition: FSQDQM.h:103
MonitorElement * PFJeteta
Definition: FSQDQM.h:107
int k[5][pyjets_maxn]
MonitorElement * CastorJetphi
Definition: FSQDQM.h:111
MonitorElement * h_leadingtrkpt_ntrk_towards
Definition: FSQDQM.h:136
unsigned int bxNumber_
Definition: FSQDQM.h:103
double x() const
x coordinate
Definition: Vertex.h:111
MonitorElement * PFJetpt
Definition: FSQDQM.h:106
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:125
double xError() const
error on x
Definition: Vertex.h:119
MonitorElement * h_leadingtrkpt_ptsum_transverse
Definition: FSQDQM.h:140
bool isFake() const
Definition: Vertex.h:72
T const * product() const
Definition: Handle.h:81
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
MonitorElement * Track_HP_Pt
Definition: FSQDQM.h:117
edm::EDGetTokenT< reco::TrackCollection > tok_track_
Definition: FSQDQM.h:98
std::vector< PFJet > PFJetCollection
collection of PFJet objects
MonitorElement * Track_HP_dzvtx_over_dzerr
Definition: FSQDQM.h:119
MonitorElement * PV_chi2
Definition: FSQDQM.h:122
fixed size matrix
HLT enums.
std::vector< TrackBaseRef >::const_iterator trackRef_iterator
The iteratator for the vector<TrackRef>
Definition: Vertex.h:37
MonitorElement * CastorJetMulti
Definition: FSQDQM.h:112
MonitorElement * h_ptsum_away
Definition: FSQDQM.h:128
double normalizedChi2() const
chi-squared divided by n.d.o.f.
Definition: Vertex.h:107
trackRef_iterator tracks_begin() const
first iterator over tracks
Definition: Vertex.cc:76
MonitorElement * h_leadingtrkpt_ntrk_away
Definition: FSQDQM.h:135
double yError() const
error on y
Definition: Vertex.h:121
size_t tracksSize() const
number of tracks
Definition: Vertex.cc:71
MonitorElement * h_ntracks_transverse
Definition: FSQDQM.h:130
MonitorElement * Track_HP_ptErr_over_pt
Definition: FSQDQM.h:118
void FSQDQM::bookHistograms ( DQMStore::IBooker bei,
edm::Run const &  ,
edm::EventSetup const &   
)
overrideprivatevirtual

Implements DQMEDAnalyzer.

Definition at line 115 of file FSQDQM.cc.

References DQMStore::IBooker::book1D(), DQMStore::IBooker::bookProfile(), and DQMStore::IBooker::setCurrentFolder().

116  {
117  bei.setCurrentFolder("Physics/FSQ");
118  PFJetpt = bei.book1D("PFJetpt",";p_{T}(PFJet)", 100,0.0 , 100);
119  PFJeteta = bei.book1D("PFJeteta", ";#eta(PFJet)", 50, -2.5, 2.5);
120  PFJetphi = bei.book1D("PFJetphi", ";#phi(PFJet)", 50, -3.14,3.14);
121  PFJetMulti = bei.book1D("PFJetMulti", ";No. of PFJets", 10, -0.5, 9.5);
122  PFJetRapidity = bei.book1D("PFJetRapidity",";PFJetRapidity", 50, -6.0,6.0);
123  CastorJetphi = bei.book1D("CastorJetphi", ";#phi(CastorJet)", 50, -3.14,3.14);
124  CastorJetMulti = bei.book1D("CastorJetMulti", ";No. of CastorJets", 10, -0.5, 9.5);
125  Track_HP_Phi =bei.book1D("Track_HP_Phi",";#phi(HPtrack)", 50, -3.14,3.14);
126  Track_HP_Eta =bei.book1D("Track_HP_Eta", ";#eta(HPtrack)", 50, -2.5, 2.5);
127  Track_HP_Pt =bei.book1D("Track_HP_Pt", ";p_{T}(HPtrack)",500, 0.0 , 50);
128  Track_HP_ptErr_over_pt=bei.book1D("Track_HP_ptErr_over_pt",";{p_{T}Err}/p_{T}",100,0,0.1);
129  Track_HP_dzvtx_over_dzerr=bei.book1D("Track_HP_dzvtx_over_dzerr",";dZerr/dZ",100,-10,10);
130  Track_HP_dxyvtx_over_dxyerror =bei.book1D("Track_HP_dxyvtx_over_dxyerror",";dxyErr/dxy",100,-10,10);
131  NPV = bei.book1D("NPV",";NPV",10, -0.5, 9.5);
132  PV_chi2 = bei.book1D("PV_chi2",";PV_chi2",100, 0.0, 2.0);
133  PV_d0 = bei.book1D("PV_d0",";PV_d0",100, -10.0, 10.0);
134  PV_numTrks = bei.book1D("PV_numTrks",";PV_numTrks",100, -0.5, 99.5);
135  PV_sumTrks=bei.book1D("PV_sumTrks",";PV_sumTrks",100,0,100);
136  h_ptsum_towards = bei.book1D("h_ptsum_towards",";h_ptsum_towards",100,0,100);
137  h_ptsum_transverse = bei.book1D("h_ptsum_transverse",";h_ptsum_transverse",100,0,100);
138 
139  h_ntracks = bei.book1D("h_ntracks",";h_ntracks",50,-0.5,49.5);
140  h_trkptsum = bei.book1D("h_trkptsum",";h_trkptsum",100,0,100);
141  h_ptsum_away = bei.book1D("h_ptsum_away",";h_ptsum_away",100,0,100);
142  h_ntracks_towards = bei.book1D("h_ntracks_towards",";h_ntracks_towards",50,-0.5,49.5);
143  h_ntracks_transverse = bei.book1D("h_ntracks_transverse",";h_ntracks_transverse",50,-0.5,49.5);
144  h_ntracks_away = bei.book1D("h_ntracks_away",";h_ntracks_away",50,-0.5,49.5);
145 
146  h_leadingtrkpt_ntrk_away =bei.bookProfile("h_leadingtrkpt_ntrk_away","h_leadingtrkpt_ntrk_away",50,0,50,0,30," ");
147  h_leadingtrkpt_ntrk_towards =bei.bookProfile("h_leadingtrkpt_ntrk_towards","h_leadingtrkpt_ntrk_towards",50,0,50,0,30," ");
148  h_leadingtrkpt_ntrk_transverse =bei.bookProfile("h_leadingtrkpt_ntrk_transverse","h_leadingtrkpt_ntrk_transverse",50,0,50,0,30," ");
149  h_leadingtrkpt_ptsum_away =bei.bookProfile("h_leadingtrkpt_ptsum_away","h_leadingtrkpt_ptsum_away",50,0,50,0,30," ");
150  h_leadingtrkpt_ptsum_towards =bei.bookProfile("h_leadingtrkpt_ptsum_towards","h_leadingtrkpt_ptsum_towards",50,0,50,0,30," ");
151  h_leadingtrkpt_ptsum_transverse =bei.bookProfile("h_leadingtrkpt_ptsum_transverse","h_leadingtrkpt_ptsum_transverse",50,0,50,0,30," ");
152 
153  }
MonitorElement * h_leadingtrkpt_ntrk_transverse
Definition: FSQDQM.h:137
MonitorElement * bookProfile(Args &&...args)
Definition: DQMStore.h:160
MonitorElement * PV_d0
Definition: FSQDQM.h:123
MonitorElement * NPV
Definition: FSQDQM.h:121
MonitorElement * PV_numTrks
Definition: FSQDQM.h:124
MonitorElement * Track_HP_dxyvtx_over_dxyerror
Definition: FSQDQM.h:120
MonitorElement * h_ptsum_transverse
Definition: FSQDQM.h:127
MonitorElement * h_ntracks_away
Definition: FSQDQM.h:131
MonitorElement * h_ptsum_towards
Definition: FSQDQM.h:126
MonitorElement * h_ntracks
Definition: FSQDQM.h:133
MonitorElement * PFJetMulti
Definition: FSQDQM.h:113
MonitorElement * h_ntracks_towards
Definition: FSQDQM.h:129
MonitorElement * Track_HP_Phi
Definition: FSQDQM.h:115
MonitorElement * PV_sumTrks
Definition: FSQDQM.h:125
MonitorElement * Track_HP_Eta
Definition: FSQDQM.h:116
MonitorElement * PFJetRapidity
Definition: FSQDQM.h:114
MonitorElement * h_trkptsum
Definition: FSQDQM.h:132
MonitorElement * h_leadingtrkpt_ptsum_away
Definition: FSQDQM.h:138
MonitorElement * PFJetphi
Definition: FSQDQM.h:108
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
MonitorElement * h_leadingtrkpt_ptsum_towards
Definition: FSQDQM.h:139
MonitorElement * PFJeteta
Definition: FSQDQM.h:107
MonitorElement * CastorJetphi
Definition: FSQDQM.h:111
MonitorElement * h_leadingtrkpt_ntrk_towards
Definition: FSQDQM.h:136
MonitorElement * PFJetpt
Definition: FSQDQM.h:106
MonitorElement * h_leadingtrkpt_ptsum_transverse
Definition: FSQDQM.h:140
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:274
MonitorElement * Track_HP_Pt
Definition: FSQDQM.h:117
MonitorElement * Track_HP_dzvtx_over_dzerr
Definition: FSQDQM.h:119
MonitorElement * PV_chi2
Definition: FSQDQM.h:122
MonitorElement * CastorJetMulti
Definition: FSQDQM.h:112
MonitorElement * h_ptsum_away
Definition: FSQDQM.h:128
MonitorElement * h_leadingtrkpt_ntrk_away
Definition: FSQDQM.h:135
MonitorElement * h_ntracks_transverse
Definition: FSQDQM.h:130
MonitorElement * Track_HP_ptErr_over_pt
Definition: FSQDQM.h:118
void FSQDQM::bookHistos ( DQMStore bei)
private

Member Data Documentation

unsigned int FSQDQM::bxNumber_
private

Definition at line 103 of file FSQDQM.h.

MonitorElement* FSQDQM::CastorJetMulti
private

Definition at line 112 of file FSQDQM.h.

MonitorElement* FSQDQM::CastorJetphi
private

Definition at line 111 of file FSQDQM.h.

unsigned int FSQDQM::eventNumber_
private

Definition at line 103 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ntrk_away
private

Definition at line 135 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ntrk_towards
private

Definition at line 136 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ntrk_transverse
private

Definition at line 137 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ptsum_away
private

Definition at line 138 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ptsum_towards
private

Definition at line 139 of file FSQDQM.h.

MonitorElement* FSQDQM::h_leadingtrkpt_ptsum_transverse
private

Definition at line 140 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ntracks
private

Definition at line 133 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ntracks_away
private

Definition at line 131 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ntracks_towards
private

Definition at line 129 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ntracks_transverse
private

Definition at line 130 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ptsum_away
private

Definition at line 128 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ptsum_towards
private

Definition at line 126 of file FSQDQM.h.

MonitorElement* FSQDQM::h_ptsum_transverse
private

Definition at line 127 of file FSQDQM.h.

MonitorElement* FSQDQM::h_trkptsum
private

Definition at line 132 of file FSQDQM.h.

std::vector<int> FSQDQM::hltresults
private

Definition at line 102 of file FSQDQM.h.

std::string FSQDQM::labelBS_
private

Definition at line 96 of file FSQDQM.h.

std::string FSQDQM::labelCastorJet_
private

Definition at line 96 of file FSQDQM.h.

std::string FSQDQM::labelPFJet_
private

Definition at line 96 of file FSQDQM.h.

std::string FSQDQM::labelTrack_
private

Definition at line 96 of file FSQDQM.h.

unsigned int FSQDQM::lumiNumber_
private

Definition at line 103 of file FSQDQM.h.

MonitorElement* FSQDQM::NPV
private

Definition at line 121 of file FSQDQM.h.

MonitorElement* FSQDQM::PFJeteta
private

Definition at line 107 of file FSQDQM.h.

MonitorElement* FSQDQM::PFJetMulti
private

Definition at line 113 of file FSQDQM.h.

MonitorElement* FSQDQM::PFJetphi
private

Definition at line 108 of file FSQDQM.h.

MonitorElement* FSQDQM::PFJetpt
private

Definition at line 106 of file FSQDQM.h.

MonitorElement* FSQDQM::PFJetRapidity
private

Definition at line 114 of file FSQDQM.h.

MonitorElement* FSQDQM::PV_chi2
private

Definition at line 122 of file FSQDQM.h.

MonitorElement* FSQDQM::PV_d0
private

Definition at line 123 of file FSQDQM.h.

MonitorElement* FSQDQM::PV_numTrks
private

Definition at line 124 of file FSQDQM.h.

MonitorElement* FSQDQM::PV_sumTrks
private

Definition at line 125 of file FSQDQM.h.

edm::EDGetTokenT<edm::View<reco::Vertex> > FSQDQM::pvs_
private

Definition at line 97 of file FSQDQM.h.

unsigned int FSQDQM::runNumber_
private

Definition at line 103 of file FSQDQM.h.

edm::EDGetTokenT<reco::BasicJetCollection> FSQDQM::tok_castorjet_
private

Definition at line 100 of file FSQDQM.h.

edm::EDGetTokenT<reco::PFJetCollection> FSQDQM::tok_pfjet_
private

Definition at line 99 of file FSQDQM.h.

edm::EDGetTokenT<reco::TrackCollection> FSQDQM::tok_track_
private

Definition at line 98 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_dxyvtx_over_dxyerror
private

Definition at line 120 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_dzvtx_over_dzerr
private

Definition at line 119 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_Eta
private

Definition at line 116 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_Phi
private

Definition at line 115 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_Pt
private

Definition at line 117 of file FSQDQM.h.

MonitorElement* FSQDQM::Track_HP_ptErr_over_pt
private

Definition at line 118 of file FSQDQM.h.

edm::InputTag FSQDQM::vertex_
private

Definition at line 95 of file FSQDQM.h.