CMS 3D CMS Logo

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

#include <DTLocalTriggerSynchTask.h>

Inheritance diagram for DTLocalTriggerSynchTask:
one::DQMEDAnalyzer< T > one::dqmimplementation::DQMBaseClass< T... >

Public Member Functions

 DTLocalTriggerSynchTask (const edm::ParameterSet &ps)
 Constructor. More...
 
 ~DTLocalTriggerSynchTask () override
 Destructor. More...
 
- Public Member Functions inherited from one::DQMEDAnalyzer< T >
 DQMEDAnalyzer ()=default
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > const &)=delete
 
 DQMEDAnalyzer (DQMEDAnalyzer< T... > &&)=delete
 
 ~DQMEDAnalyzer () override=default
 

Protected Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &context) override
 Analyze. More...
 
std::string & baseDir ()
 
void bookHistograms (DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
 Book the histograms. More...
 
void bookHistos (DQMStore::IBooker &, const DTChamberId &dtCh)
 Book the histograms. More...
 
void dqmBeginRun (const edm::Run &, const edm::EventSetup &) override
 Beginrun. More...
 

Protected Attributes

const int wheelArrayShift = 3
 

Private Attributes

float angleRange
 
std::string baseDirectory
 
float bxTime
 
edm::EDGetTokenT< DTLocalTriggerCollectionddu_Token_
 
int fineDelay
 
float minHitsPhi
 
edm::ESHandle< DTGeometrymuonGeom
 
int nBXHigh
 
int nBXLow
 
int nevents
 
DTArr3int phCodeBestDDU
 
DTArr3int phCodeBestTM
 
DTArr4int phCodeBXTM
 
bool processDDU
 
bool rangeInBX
 
edm::EDGetTokenT< DTRecSegment4DCollectionseg_Token_
 
DTArr3int segHitBest
 
DTArr3int thCodeBestDDU
 
MonitorElementtm_IDDataErrorPlot
 
edm::EDGetTokenT< L1MuDTChambPhContainertm_Token_
 
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
 
std::unique_ptr< DTTTrigBaseSynctTrigSync
 

Friends

class DTMonitorModule
 

Detailed Description

Definition at line 51 of file DTLocalTriggerSynchTask.h.

Constructor & Destructor Documentation

DTLocalTriggerSynchTask::DTLocalTriggerSynchTask ( const edm::ParameterSet ps)

Constructor.

Definition at line 43 of file DTLocalTriggerSynchTask.cc.

References angleRange, baseDirectory, bxTime, ddu_Token_, reco::get(), edm::ParameterSet::getParameter(), minHitsPhi, nBXHigh, nBXLow, processDDU, rangeInBX, seg_Token_, AlCaHLTBitMon_QueryRunRegistry::string, and tm_Token_.

43  : nevents(0),
44  tTrigSync{DTTTrigSyncFactory::get()->create(ps.getParameter<std::string>("tTrigMode"),
45  ps.getParameter<edm::ParameterSet>("tTrigModeConfig"))}
46 {
47 
48 
49  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Constructor" << endl;
50  tm_Token_ = consumes<L1MuDTChambPhContainer>(
51  ps.getParameter<edm::InputTag>("TMInputTag"));
52  seg_Token_ = consumes<DTRecSegment4DCollection>(
53  ps.getParameter<edm::InputTag>("SEGInputTag"));
54 
55  processDDU = ps.getUntrackedParameter<bool>("processDDU",false);
56 
57  if (processDDU)
58  ddu_Token_ = consumes<DTLocalTriggerCollection>(
59  ps.getParameter<edm::InputTag>("DDUInputTag"));
60 
61  bxTime = ps.getParameter<double>("bxTimeInterval"); // CB move this to static const or DB
62  rangeInBX = ps.getParameter<bool>("rangeWithinBX");
63  nBXLow = ps.getParameter<int>("nBXLow");
64  nBXHigh = ps.getParameter<int>("nBXHigh");
65  angleRange = ps.getParameter<double>("angleRange");
66  minHitsPhi = ps.getParameter<int>("minHitsPhi");
67  baseDirectory = ps.getParameter<string>("baseDir");
68 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< DTTTrigBaseSync > tTrigSync
edm::EDGetTokenT< DTLocalTriggerCollection > ddu_Token_
edm::EDGetTokenT< DTRecSegment4DCollection > seg_Token_
edm::EDGetTokenT< L1MuDTChambPhContainer > tm_Token_
T get(const Candidate &c)
Definition: component.h:55
DTLocalTriggerSynchTask::~DTLocalTriggerSynchTask ( )
override

Destructor.

Definition at line 71 of file DTLocalTriggerSynchTask.cc.

References nevents.

71  {
72 
73  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
74 
75 }

Member Function Documentation

void DTLocalTriggerSynchTask::analyze ( const edm::Event event,
const edm::EventSetup context 
)
overrideprotected

Analyze.

Definition at line 106 of file DTLocalTriggerSynchTask.cc.

References angleRange, bxTime, ddu_Token_, DTRecSegment2D::degreesOfFreedom(), dir, L1MuDTChambPhContainer::getContainer(), hfClusterShapes_cfi::hits, mps_fire::i, createfilelist::int, DTRecSegment2D::ist0Valid(), gen::k, minHitsPhi, nevents, phCodeBestDDU, phCodeBestTM, Geom::pi(), processDDU, jets_cff::quality, rangeInBX, DetId::rawId(), DTChamberId::sector(), seg_Token_, DTChamberId::station(), relativeConstraints::station, ntuplemaker::time, tm_Token_, HiIsolationCommonParameters_cff::track, triggerHistos, tTrigSync, DTChamberId::wheel(), makeMuonMisalignmentScenario::wheel, and wheelArrayShift.

106  {
107 
108  nevents++;
109 
110  for (int i=0;i<5;++i){
111  for (int j=0;j<6;++j){
112  for (int k=0;k<13;++k){
113  phCodeBestTM[j][i][k] = -1;
114  phCodeBestDDU[j][i][k] = -1;
115  }
116  }
117  }
118 
119  // Get best TM triggers
121  event.getByToken(tm_Token_, l1DTTPGPh);
122  vector<L1MuDTChambPhDigi> const* phTrigs = l1DTTPGPh->getContainer();
123 
124  vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin();
125  vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
126  for(; iph !=iphe ; ++iph) {
127 
128  int phwheel = iph->whNum();
129  int phsec = iph->scNum() + 1; // DTTF[0-11] -> DT[1-12] Sector Numbering
130  int phst = iph->stNum();
131  int phcode = iph->code();
132 
133  if(phcode>phCodeBestTM[phwheel+3][phst][phsec] && phcode<7) {
134  phCodeBestTM[phwheel+3][phst][phsec]=phcode;
135  }
136 
137  }
138 
139  // Get best DDU triggers
140  if (processDDU){
142  event.getByToken(ddu_Token_, trigsDDU);
144 
145  for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){
146 
147  const DTChamberId& id = (*detUnitIt).first;
148  const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
149 
150  int wh = id.wheel();
151  int sec = id.sector();
152  int st = id.station();
153 
154  for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt!=range.second;++trigIt){
155 
156  int quality = trigIt->quality();
157 
158  if(quality>-1 && quality<7 &&
159  quality>phCodeBestDDU[wh+wheelArrayShift][st][sec]) {
161  }
162  }
163  }
164  }
165 
166  //Get best segments (highest number of phi hits)
167  vector<const DTRecSegment4D*> bestSegments4D;
169  event.getByToken(seg_Token_, segments4D);
172 
173  for (chambIdIt = segments4D->id_begin(); chambIdIt != segments4D->id_end(); ++chambIdIt){
174 
175  DTRecSegment4DCollection::range range = segments4D->get(*chambIdIt);
176  const DTRecSegment4D* best=nullptr;
177  int hitsBest = 0;
178  int hits = 0;
179 
180  for ( track = range.first; track != range.second; ++track){
181  if( (*track).hasPhi() ) {
182  hits = (*track).phiSegment()->degreesOfFreedom()+2;
183  if ( hits>hitsBest ){
184  best = &(*track);
185  hitsBest = hits;
186  }
187  }
188  }
189  if (best) {
190  bestSegments4D.push_back(best);
191  }
192  }
193 
194 
195  // Filling histos
196  vector<const DTRecSegment4D*>::const_iterator bestSegIt = bestSegments4D.begin();
197  vector<const DTRecSegment4D*>::const_iterator bestSegEnd = bestSegments4D.end();
198  for (; bestSegIt!=bestSegEnd; ++bestSegIt ){
199 
200  float dir = atan((*bestSegIt)->localDirection().x()/ (*bestSegIt)->localDirection().z())*180/Geom::pi(); // CB cerca un modo migliore x farlo
201  const DTRecSegment2D* seg2D = (*bestSegIt)->phiSegment();
202  int nHitsPhi = seg2D->degreesOfFreedom()+2;
203  DTChamberId chambId = (*bestSegIt)->chamberId();
204  map<string, MonitorElement*> &innerME = triggerHistos[chambId.rawId()];
205 
206  if (fabs(dir)<angleRange &&
207  nHitsPhi>=minHitsPhi &&
208  seg2D->ist0Valid()){
209 
210  float t0seg = (*bestSegIt)->phiSegment()->t0();
211  float tTrig = (tTrigSync->offset(DTWireId(chambId,1,1,2)) + tTrigSync->offset(DTWireId(chambId,3,1,2)) )/2;
212  float time = tTrig+t0seg;
213  float htime = rangeInBX ? time-int(time/bxTime)*bxTime : time-int(tTrig/bxTime)*bxTime;
214 
215  int wheel = chambId.wheel();
216  int sector = chambId.sector();
217  int station = chambId.station();
218  int scsector = sector>12 ? sector==13 ? 4 : 10 : sector;
219 
220  int qualTM = phCodeBestTM[wheel+3][station][scsector];
221  int qualDDU = phCodeBestDDU[wheel+3][station][scsector];
222 
223  if (fabs(t0seg)>0.01) {
224  innerME.find("SEG_TrackCrossingTime")->second->Fill(htime);
225  if ( qualTM>=0 ) innerME.find("TM_TrackCrossingTimeAll")->second->Fill(htime);
226  if ( qualTM==6 ) innerME.find("TM_TrackCrossingTimeHH")->second->Fill(htime);
227  if ( processDDU && qualDDU>=0 ) innerME.find("DDU_TrackCrossingTimeAll")->second->Fill(htime);
228  if ( processDDU && qualDDU==6 ) innerME.find("DDU_TrackCrossingTimeHH")->second->Fill(htime);
229  }
230 
231  }
232  }
233 
234 }
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:50
identifier iterator
Definition: RangeMap.h:135
std::unique_ptr< DTTTrigBaseSync > tTrigSync
edm::EDGetTokenT< DTLocalTriggerCollection > ddu_Token_
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
int degreesOfFreedom() const override
return the DOF of the segment
bool ist0Valid() const
edm::EDGetTokenT< DTRecSegment4DCollection > seg_Token_
int k[5][pyjets_maxn]
std::vector< DigiType >::const_iterator const_iterator
Phi_Container const * getContainer() const
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
int sector() const
Definition: DTChamberId.h:61
std::pair< const_iterator, const_iterator > Range
dbl *** dir
Definition: mlp_gen.cc:35
constexpr double pi()
Definition: Pi.h:31
int station() const
Return the station number.
Definition: DTChamberId.h:51
edm::EDGetTokenT< L1MuDTChambPhContainer > tm_Token_
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
std::string& DTLocalTriggerSynchTask::baseDir ( )
inlineprotected

Definition at line 77 of file DTLocalTriggerSynchTask.h.

References baseDirectory.

Referenced by bookHistograms(), and bookHistos().

77 { return baseDirectory; }
void DTLocalTriggerSynchTask::bookHistograms ( DQMStore::IBooker ibooker,
edm::Run const &  iRun,
edm::EventSetup const &  context 
)
overrideprotected

Book the histograms.

Definition at line 77 of file DTLocalTriggerSynchTask.cc.

References baseDir(), DQMStore::IBooker::bookFloat(), bookHistos(), bxTime, DTGeometry::chambers(), MonitorElement::Fill(), muonGeom, DQMStore::IBooker::setCurrentFolder(), triggerHistos, and tTrigSync.

77  {
78 
79  edm::LogVerbatim ("DTLocalTriggerSynchTask") <<"[DTLocalTriggerSynchTask]: Book Histograms"<<endl;
80 
81  ibooker.setCurrentFolder(baseDir());
82  ibooker.bookFloat("BXTimeSpacing")->Fill(bxTime);
83 
84  tTrigSync->setES(context);
85 
86 
87  std::vector<const DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin();
88  std::vector<const DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();
89 
90  for (; chambIt!=chambEnd; ++chambIt) {
91  bookHistos(ibooker,(*chambIt)->id());
92  triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL1"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),1,1,2)));
93  triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL3"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),3,1,2)));
94  }
95 
96 
97 }
const std::vector< const DTChamber * > & chambers() const
Return a vector of all Chamber.
Definition: DTGeometry.cc:102
std::unique_ptr< DTTTrigBaseSync > tTrigSync
void Fill(long long x)
edm::ESHandle< DTGeometry > muonGeom
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
void bookHistos(DQMStore::IBooker &, const DTChamberId &dtCh)
Book the histograms.
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:105
void DTLocalTriggerSynchTask::bookHistos ( DQMStore::IBooker ibooker,
const DTChamberId dtCh 
)
protected

Book the histograms.

Definition at line 236 of file DTLocalTriggerSynchTask.cc.

References baseDir(), DQMStore::IBooker::book1D(), DQMStore::IBooker::bookFloat(), bxTime, SiStripPI::max, min(), pileupCalc::nbins, nBXHigh, nBXLow, processDDU, rangeInBX, DetId::rawId(), DTChamberId::sector(), DQMStore::IBooker::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, triggerHistos, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by bookHistograms().

236  {
237 
238  stringstream wheel; wheel << dtChId.wheel();
239  stringstream station; station << dtChId.station();
240  stringstream sector; sector << dtChId.sector();
241  uint32_t chRawId = dtChId.rawId();
242 
243  ibooker.setCurrentFolder(baseDir() + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() );
244 
245  std::vector<string> histoTags = { "SEG_TrackCrossingTime", "TM_TrackCrossingTimeAll", "TM_TrackCrossingTimeHH" };
246 
247  if (processDDU) {
248  histoTags.push_back("DDU_TrackCrossingTimeAll");
249  histoTags.push_back("DDU_TrackCrossingTimeHH");
250  }
251 
252  float min = rangeInBX ? 0 : nBXLow*bxTime;
253  float max = rangeInBX ? bxTime : nBXHigh*bxTime;
254  int nbins = static_cast<int>(ceil( rangeInBX ? bxTime : (nBXHigh-nBXLow)*bxTime));
255 
256  for (const auto & histoTag : histoTags) {
257  string histoName = histoTag + (rangeInBX ? "InBX" : "") + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
258  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: booking "
259  << baseDir() + "/Wheel" << wheel.str()
260  << "/Sector" << sector.str()
261  << "/Station"<< station.str()
262  << "/" << histoName << endl;
263 
264  triggerHistos[chRawId][histoTag] = ibooker.book1D(histoName.c_str(),"Track time distribution",nbins,min,max);
265  }
266 
267  string floatTag[2] = { "tTrig_SL1", "tTrig_SL3" };
268 
269  for (int iFloat=0;iFloat<2;++iFloat) {
270  string floatName = floatTag[iFloat] + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
271  triggerHistos[chRawId][floatTag[iFloat]] = ibooker.bookFloat(floatName);
272  }
273 
274 }
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
T min(T a, T b)
Definition: MathUtil.h:58
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
MonitorElement * bookFloat(Args &&...args)
Definition: DQMStore.h:105
void DTLocalTriggerSynchTask::dqmBeginRun ( const edm::Run run,
const edm::EventSetup context 
)
overrideprotected

Beginrun.

Definition at line 100 of file DTLocalTriggerSynchTask.cc.

References edm::EventSetup::get(), and muonGeom.

100  {
101 
102  context.get<MuonGeometryRecord>().get(muonGeom);
103 }
edm::ESHandle< DTGeometry > muonGeom
T get() const
Definition: EventSetup.h:71

Friends And Related Function Documentation

friend class DTMonitorModule
friend

Definition at line 53 of file DTLocalTriggerSynchTask.h.

Member Data Documentation

float DTLocalTriggerSynchTask::angleRange
private

Definition at line 96 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and DTLocalTriggerSynchTask().

std::string DTLocalTriggerSynchTask::baseDirectory
private

Definition at line 101 of file DTLocalTriggerSynchTask.h.

Referenced by baseDir(), and DTLocalTriggerSynchTask().

float DTLocalTriggerSynchTask::bxTime
private
edm::EDGetTokenT<DTLocalTriggerCollection> DTLocalTriggerSynchTask::ddu_Token_
private

Definition at line 108 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and DTLocalTriggerSynchTask().

int DTLocalTriggerSynchTask::fineDelay
private

Definition at line 98 of file DTLocalTriggerSynchTask.h.

float DTLocalTriggerSynchTask::minHitsPhi
private

Definition at line 97 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and DTLocalTriggerSynchTask().

edm::ESHandle<DTGeometry> DTLocalTriggerSynchTask::muonGeom
private

Definition at line 103 of file DTLocalTriggerSynchTask.h.

Referenced by bookHistograms(), and dqmBeginRun().

int DTLocalTriggerSynchTask::nBXHigh
private

Definition at line 95 of file DTLocalTriggerSynchTask.h.

Referenced by bookHistos(), and DTLocalTriggerSynchTask().

int DTLocalTriggerSynchTask::nBXLow
private

Definition at line 94 of file DTLocalTriggerSynchTask.h.

Referenced by bookHistos(), and DTLocalTriggerSynchTask().

int DTLocalTriggerSynchTask::nevents
private

Definition at line 83 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and ~DTLocalTriggerSynchTask().

DTArr3int DTLocalTriggerSynchTask::phCodeBestDDU
private

Definition at line 87 of file DTLocalTriggerSynchTask.h.

Referenced by analyze().

DTArr3int DTLocalTriggerSynchTask::phCodeBestTM
private

Definition at line 85 of file DTLocalTriggerSynchTask.h.

Referenced by analyze().

DTArr4int DTLocalTriggerSynchTask::phCodeBXTM
private

Definition at line 86 of file DTLocalTriggerSynchTask.h.

bool DTLocalTriggerSynchTask::processDDU
private

Definition at line 93 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), bookHistos(), and DTLocalTriggerSynchTask().

bool DTLocalTriggerSynchTask::rangeInBX
private

Definition at line 92 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), bookHistos(), and DTLocalTriggerSynchTask().

edm::EDGetTokenT<DTRecSegment4DCollection> DTLocalTriggerSynchTask::seg_Token_
private

Definition at line 109 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and DTLocalTriggerSynchTask().

DTArr3int DTLocalTriggerSynchTask::segHitBest
private

Definition at line 89 of file DTLocalTriggerSynchTask.h.

DTArr3int DTLocalTriggerSynchTask::thCodeBestDDU
private

Definition at line 88 of file DTLocalTriggerSynchTask.h.

MonitorElement* DTLocalTriggerSynchTask::tm_IDDataErrorPlot
private

Definition at line 105 of file DTLocalTriggerSynchTask.h.

edm::EDGetTokenT<L1MuDTChambPhContainer> DTLocalTriggerSynchTask::tm_Token_
private

Definition at line 107 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and DTLocalTriggerSynchTask().

std::map<uint32_t, std::map<std::string, MonitorElement*> > DTLocalTriggerSynchTask::triggerHistos
private

Definition at line 104 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), bookHistograms(), and bookHistos().

std::unique_ptr<DTTTrigBaseSync> DTLocalTriggerSynchTask::tTrigSync
private

Definition at line 99 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and bookHistograms().

const int DTLocalTriggerSynchTask::wheelArrayShift = 3
protected

Definition at line 79 of file DTLocalTriggerSynchTask.h.

Referenced by analyze().