CMS 3D CMS Logo

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

#include <DTLocalTriggerSynchTask.h>

Inheritance diagram for DTLocalTriggerSynchTask:
edm::EDAnalyzer

Public Member Functions

 DTLocalTriggerSynchTask (const edm::ParameterSet &ps)
 Constructor. More...
 
virtual ~DTLocalTriggerSynchTask ()
 Destructor. More...
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Protected Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &context)
 Analyze. More...
 
std::string & baseDir ()
 
void beginJob ()
 
void beginRun (const edm::Run &run, const edm::EventSetup &context)
 Begin Run. More...
 
void bookHistos (const DTChamberId &dtCh)
 Book the histograms. More...
 
void endJob (void)
 EndJob. More...
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Private Attributes

float angleRange
 
std::string baseDirectory
 
float bxTime
 
DQMStoredbe
 
MonitorElementdcc_IDDataErrorPlot
 
int fineDelay
 
float minHitsPhi
 
edm::ESHandle< DTGeometrymuonGeom
 
int nBXHigh
 
int nBXLow
 
int nevents
 
edm::ParameterSet parameters
 
int phCodeBestDCC [6][5][13]
 
int phCodeBestDDU [6][5][13]
 
int phCodeBXDCC [6][5][13][3]
 
bool rangeInBX
 
int segHitBest [6][5][13]
 
int thCodeBestDDU [6][5][13]
 
std::map< uint32_t, std::map
< std::string, MonitorElement * > > 
triggerHistos
 
DTTTrigBaseSynctTrigSync
 

Friends

class DTMonitorModule
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
typedef WorkerT< EDAnalyzerWorkerType
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 

Detailed Description

Definition at line 44 of file DTLocalTriggerSynchTask.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 50 of file DTLocalTriggerSynchTask.cc.

References parameters.

50  : nevents(0) {
51 
52  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: Constructor" << endl;
53  parameters = ps;
54 
55 }
DTLocalTriggerSynchTask::~DTLocalTriggerSynchTask ( )
virtual

Destructor.

Definition at line 58 of file DTLocalTriggerSynchTask.cc.

References nevents.

58  {
59 
60  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
61 
62 }

Member Function Documentation

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

Analyze.

Implements edm::EDAnalyzer.

Definition at line 113 of file DTLocalTriggerSynchTask.cc.

References angleRange, bxTime, DTRecSegment2D::degreesOfFreedom(), dir, edm::ParameterSet::getParameter(), i, DTRecSegment2D::ist0Valid(), j, gen::k, minHitsPhi, nevents, DTTTrigBaseSync::offset(), parameters, phCodeBestDCC, phCodeBestDDU, DTRecSegment4D::phiSegment(), Geom::pi(), rangeInBX, DetId::rawId(), DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, cond::rpcobgas::time, triggerHistos, tTrigSync, and DTChamberId::wheel().

113  {
114 
115  nevents++;
116 
117  InputTag inputTagDCC = parameters.getParameter<edm::InputTag>("DCCInputTag");
118  InputTag inputTagDDU = parameters.getParameter<edm::InputTag>("DDUInputTag");
119  InputTag inputTagSEG = parameters.getParameter<edm::InputTag>("SEGInputTag");
120 
121  for (int i=0;i<5;++i){
122  for (int j=0;j<6;++j){
123  for (int k=0;k<13;++k){
124  phCodeBestDCC[j][i][k] = -1;
125  phCodeBestDDU[j][i][k] = -1;
126  }
127  }
128  }
129 
130  // Get best DCC triggers
132  event.getByLabel(inputTagDCC,l1DTTPGPh);
133  vector<L1MuDTChambPhDigi>* phTrigs = l1DTTPGPh->getContainer();
134 
135  vector<L1MuDTChambPhDigi>::const_iterator iph = phTrigs->begin();
136  vector<L1MuDTChambPhDigi>::const_iterator iphe = phTrigs->end();
137  for(; iph !=iphe ; ++iph) {
138 
139  int phwheel = iph->whNum();
140  int phsec = iph->scNum() + 1; // DTTF[0-11] -> DT[1-12] Sector Numbering
141  int phst = iph->stNum();
142  int phcode = iph->code();
143 
144  if(phcode>phCodeBestDCC[phwheel+3][phst][phsec] && phcode<7) {
145  phCodeBestDCC[phwheel+3][phst][phsec]=phcode;
146  }
147 
148  }
149 
150  // Get best DDU triggers
152  event.getByLabel(inputTagDDU,trigsDDU);
154 
155  for (detUnitIt=trigsDDU->begin();detUnitIt!=trigsDDU->end();++detUnitIt){
156 
157  const DTChamberId& id = (*detUnitIt).first;
158  const DTLocalTriggerCollection::Range& range = (*detUnitIt).second;
159 
160  int wh = id.wheel();
161  int sec = id.sector();
162  int st = id.station();
163 
164  for (DTLocalTriggerCollection::const_iterator trigIt = range.first; trigIt!=range.second;++trigIt){
165 
166  int quality = trigIt->quality();
167 
168  if(quality>-1 && quality<7 &&
169  quality>phCodeBestDDU[wh+3][st][sec]) {
170  phCodeBestDDU[wh+3][st][sec]=quality;
171  }
172  }
173  }
174 
175  //Get best segments (highest number of phi hits)
176  vector<const DTRecSegment4D*> bestSegments4D;
178  event.getByLabel(inputTagSEG, segments4D);
181 
182  for (chambIdIt = segments4D->id_begin(); chambIdIt != segments4D->id_end(); ++chambIdIt){
183 
184  DTRecSegment4DCollection::range range = segments4D->get(*chambIdIt);
185  const DTRecSegment4D* best=0;
186  int hitsBest = 0;
187  int hits = 0;
188 
189  for ( track = range.first; track != range.second; ++track){
190  if( (*track).hasPhi() ) {
191  hits = (*track).phiSegment()->degreesOfFreedom()+2;
192  if ( hits>hitsBest ){
193  best = &(*track);
194  hitsBest = hits;
195  }
196  }
197  }
198  if (best) {
199  bestSegments4D.push_back(best);
200  }
201  }
202 
203 
204  // Filling histos
205  vector<const DTRecSegment4D*>::const_iterator bestSegIt = bestSegments4D.begin();
206  vector<const DTRecSegment4D*>::const_iterator bestSegEnd = bestSegments4D.end();
207  for (; bestSegIt!=bestSegEnd; ++bestSegIt ){
208 
209  float dir = atan((*bestSegIt)->localDirection().x()/ (*bestSegIt)->localDirection().z())*180/Geom::pi(); // CB cerca un modo migliore x farlo
210  const DTRecSegment2D* seg2D = (*bestSegIt)->phiSegment();
211  int nHitsPhi = seg2D->degreesOfFreedom()+2;
212  DTChamberId chambId = (*bestSegIt)->chamberId();
213  map<string, MonitorElement*> &innerME = triggerHistos[chambId.rawId()];
214 
215  if (fabs(dir)<angleRange &&
216  nHitsPhi>=minHitsPhi &&
217  seg2D->ist0Valid()){
218 
219  float t0seg = (*bestSegIt)->phiSegment()->t0();
220  float tTrig = (tTrigSync->offset(DTWireId(chambId,1,1,2)) + tTrigSync->offset(DTWireId(chambId,3,1,2)) )/2;
221  float time = tTrig+t0seg;
222  float htime = rangeInBX ? time-int(time/bxTime)*bxTime : time-int(tTrig/bxTime)*bxTime;
223 
224  int wheel = chambId.wheel();
225  int sector = chambId.sector();
226  int station = chambId.station();
227  int scsector = sector>12 ? sector==13 ? 4 : 10 : sector;
228 
229  int qualDCC = phCodeBestDCC[wheel+3][station][scsector];
230  int qualDDU = phCodeBestDDU[wheel+3][station][scsector];
231 
232  if (fabs(t0seg)>0.01) {
233  innerME.find("SEG_TrackCrossingTime")->second->Fill(htime);
234  if ( qualDCC>=0 ) innerME.find("DCC_TrackCrossingTimeAll")->second->Fill(htime);
235  if ( qualDCC==6 ) innerME.find("DCC_TrackCrossingTimeHH")->second->Fill(htime);
236  if ( qualDDU>=0 ) innerME.find("DDU_TrackCrossingTimeAll")->second->Fill(htime);
237  if ( qualDDU==6 ) innerME.find("DDU_TrackCrossingTimeHH")->second->Fill(htime);
238  }
239 
240  }
241  }
242 
243 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:53
virtual int degreesOfFreedom() const
return the DOF of the segment
const DTChamberRecSegment2D * phiSegment() const
The superPhi segment: 0 if no phi projection available.
double offset(const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globalPos)
identifier iterator
Definition: RangeMap.h:139
uint32_t rawId() const
get the raw id
Definition: DetId.h:45
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:46
int j
Definition: DBlmapReader.cc:9
bool ist0Valid() const
int k[5][pyjets_maxn]
std::vector< DigiType >::const_iterator const_iterator
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
double pi()
Definition: Pi.h:31
int sector() const
Definition: DTChamberId.h:63
std::pair< const_iterator, const_iterator > Range
dbl *** dir
Definition: mlp_gen.cc:35
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
std::string& DTLocalTriggerSynchTask::baseDir ( )
inlineprotected

Definition at line 73 of file DTLocalTriggerSynchTask.h.

References baseDirectory.

Referenced by beginJob(), bookHistos(), and endJob().

73 { return baseDirectory; }
void DTLocalTriggerSynchTask::beginJob ( void  )
protectedvirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 65 of file DTLocalTriggerSynchTask.cc.

References angleRange, baseDir(), baseDirectory, bxTime, dbe, edm::ParameterSet::getParameter(), minHitsPhi, nBXHigh, nBXLow, cmsCodeRules.cppFunctionSkipper::operator, parameters, and rangeInBX.

65  {
66 
67  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: BeginJob" << endl;
68 
69  bxTime = parameters.getParameter<double>("bxTimeInterval"); // CB move this to static const or DB
70  rangeInBX = parameters.getParameter<bool>("rangeWithinBX");
71  nBXLow = parameters.getParameter<int>("nBXLow");
72  nBXHigh = parameters.getParameter<int>("nBXHigh");
73  angleRange = parameters.getParameter<double>("angleRange");
74  minHitsPhi = parameters.getParameter<int>("minHitsPhi");
75  baseDirectory = parameters.getParameter<string>("baseDir");
76 
79  dbe->bookFloat("BXTimeSpacing")->Fill(bxTime);
80 
81 }
T getParameter(std::string const &) const
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:451
void Fill(long long x)
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void DTLocalTriggerSynchTask::beginRun ( const edm::Run run,
const edm::EventSetup context 
)
protectedvirtual

Begin Run.

Reimplemented from edm::EDAnalyzer.

Definition at line 83 of file DTLocalTriggerSynchTask.cc.

References bookHistos(), edm::EventSetup::get(), reco::get(), edm::ParameterSet::getParameter(), muonGeom, DTTTrigBaseSync::offset(), parameters, DTTTrigBaseSync::setES(), triggerHistos, and tTrigSync.

83  {
84 
85  edm::LogVerbatim ("DTLocalTriggerSynchTask") <<"[DTLocalTriggerSynchTask]: Begin Run"<<endl;
86 
87  context.get<MuonGeometryRecord>().get(muonGeom);
88  tTrigSync = DTTTrigSyncFactory::get()->create(parameters.getParameter<std::string>("tTrigMode"),
89  parameters.getParameter<edm::ParameterSet>("tTrigModeConfig"));
90  tTrigSync->setES(context);
91 
92 
93  std::vector<DTChamber*>::const_iterator chambIt = muonGeom->chambers().begin();
94  std::vector<DTChamber*>::const_iterator chambEnd = muonGeom->chambers().end();
95 
96  for (; chambIt!=chambEnd; ++chambIt) {
97  bookHistos((*chambIt)->id());
98  triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL1"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),1,1,2)));
99  triggerHistos[(*chambIt)->id().rawId()]["tTrig_SL3"]->Fill(tTrigSync->offset(DTWireId((*chambIt)->id(),3,1,2)));
100  }
101 
102 }
T getParameter(std::string const &) const
double offset(const DTLayer *layer, const DTWireId &wireId, const GlobalPoint &globalPos)
edm::ESHandle< DTGeometry > muonGeom
virtual void setES(const edm::EventSetup &setup)=0
Pass the Event Setup to the synchronization module at each event.
const T & get() const
Definition: EventSetup.h:55
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
T get(const Candidate &c)
Definition: component.h:56
void bookHistos(const DTChamberId &dtCh)
Book the histograms.
void DTLocalTriggerSynchTask::bookHistos ( const DTChamberId dtCh)
protected

Book the histograms.

Definition at line 245 of file DTLocalTriggerSynchTask.cc.

References baseDir(), DQMStore::book1D(), DQMStore::bookFloat(), bxTime, dbe, max(), min, RecoTauCommonJetSelections_cfi::nbins, nBXHigh, nBXLow, rangeInBX, DetId::rawId(), DTChamberId::sector(), DQMStore::setCurrentFolder(), DTChamberId::station(), relativeConstraints::station, triggerHistos, and DTChamberId::wheel().

Referenced by beginRun().

245  {
246 
247  stringstream wheel; wheel << dtChId.wheel();
248  stringstream station; station << dtChId.station();
249  stringstream sector; sector << dtChId.sector();
250  uint32_t chRawId = dtChId.rawId();
251 
252  dbe->setCurrentFolder(baseDir() + "/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str() );
253 
254 
255  string histoTag[5] = { "SEG_TrackCrossingTime", "DCC_TrackCrossingTimeAll", "DCC_TrackCrossingTimeHH", "DDU_TrackCrossingTimeAll", "DDU_TrackCrossingTimeHH" };
256 
257  float min = rangeInBX ? 0 : nBXLow*bxTime;
258  float max = rangeInBX ? bxTime : nBXHigh*bxTime;
259  int nbins = static_cast<int>(ceil( rangeInBX ? bxTime : (nBXHigh-nBXLow)*bxTime));
260 
261  for (int iHisto=0;iHisto<5;++iHisto) {
262  string histoName = histoTag[iHisto] + (rangeInBX ? "InBX" : "") + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
263  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: booking "
264  << baseDir() + "/Wheel" << wheel.str()
265  << "/Sector" << sector.str()
266  << "/Station"<< station.str()
267  << "/" << histoName << endl;
268 
269  triggerHistos[chRawId][histoTag[iHisto]] = dbe->book1D(histoName.c_str(),"Track time distribution",nbins,min,max);
270  }
271 
272  string floatTag[2] = { "tTrig_SL1", "tTrig_SL3" };
273 
274  for (int iFloat=0;iFloat<2;++iFloat) {
275  string floatName = floatTag[iFloat] + "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
276  triggerHistos[chRawId][floatTag[iFloat]] = dbe->bookFloat(floatName);
277  }
278 
279 }
MonitorElement * book1D(const char *name, const char *title, int nchX, double lowX, double highX)
Book 1D histogram.
Definition: DQMStore.cc:514
#define min(a, b)
Definition: mlp_lapack.h:161
MonitorElement * bookFloat(const char *name)
Book float.
Definition: DQMStore.cc:451
const T & max(const T &a, const T &b)
std::map< uint32_t, std::map< std::string, MonitorElement * > > triggerHistos
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:232
void DTLocalTriggerSynchTask::endJob ( void  )
protectedvirtual

EndJob.

Reimplemented from edm::EDAnalyzer.

Definition at line 105 of file DTLocalTriggerSynchTask.cc.

References baseDir(), dbe, nevents, and DQMStore::rmdir().

105  {
106 
107  edm::LogVerbatim ("DTLocalTriggerSynchTask") << "[DTLocalTriggerSynchTask]: analyzed " << nevents << " events" << endl;
108  dbe->rmdir(baseDir());
109 
110 }
void rmdir(const std::string &fullpath)
Definition: DQMStore.cc:2311

Friends And Related Function Documentation

friend class DTMonitorModule
friend

Definition at line 46 of file DTLocalTriggerSynchTask.h.

Member Data Documentation

float DTLocalTriggerSynchTask::angleRange
private

Definition at line 89 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and beginJob().

std::string DTLocalTriggerSynchTask::baseDirectory
private

Definition at line 94 of file DTLocalTriggerSynchTask.h.

Referenced by baseDir(), and beginJob().

float DTLocalTriggerSynchTask::bxTime
private

Definition at line 85 of file DTLocalTriggerSynchTask.h.

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

DQMStore* DTLocalTriggerSynchTask::dbe
private

Definition at line 96 of file DTLocalTriggerSynchTask.h.

Referenced by beginJob(), bookHistos(), and endJob().

MonitorElement* DTLocalTriggerSynchTask::dcc_IDDataErrorPlot
private

Definition at line 100 of file DTLocalTriggerSynchTask.h.

int DTLocalTriggerSynchTask::fineDelay
private

Definition at line 91 of file DTLocalTriggerSynchTask.h.

float DTLocalTriggerSynchTask::minHitsPhi
private

Definition at line 90 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and beginJob().

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

Definition at line 98 of file DTLocalTriggerSynchTask.h.

Referenced by beginRun().

int DTLocalTriggerSynchTask::nBXHigh
private

Definition at line 88 of file DTLocalTriggerSynchTask.h.

Referenced by beginJob(), and bookHistos().

int DTLocalTriggerSynchTask::nBXLow
private

Definition at line 87 of file DTLocalTriggerSynchTask.h.

Referenced by beginJob(), and bookHistos().

int DTLocalTriggerSynchTask::nevents
private

Definition at line 77 of file DTLocalTriggerSynchTask.h.

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

edm::ParameterSet DTLocalTriggerSynchTask::parameters
private
int DTLocalTriggerSynchTask::phCodeBestDCC[6][5][13]
private

Definition at line 79 of file DTLocalTriggerSynchTask.h.

Referenced by analyze().

int DTLocalTriggerSynchTask::phCodeBestDDU[6][5][13]
private

Definition at line 81 of file DTLocalTriggerSynchTask.h.

Referenced by analyze().

int DTLocalTriggerSynchTask::phCodeBXDCC[6][5][13][3]
private

Definition at line 80 of file DTLocalTriggerSynchTask.h.

bool DTLocalTriggerSynchTask::rangeInBX
private

Definition at line 86 of file DTLocalTriggerSynchTask.h.

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

int DTLocalTriggerSynchTask::segHitBest[6][5][13]
private

Definition at line 83 of file DTLocalTriggerSynchTask.h.

int DTLocalTriggerSynchTask::thCodeBestDDU[6][5][13]
private

Definition at line 82 of file DTLocalTriggerSynchTask.h.

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

Definition at line 99 of file DTLocalTriggerSynchTask.h.

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

DTTTrigBaseSync* DTLocalTriggerSynchTask::tTrigSync
private

Definition at line 92 of file DTLocalTriggerSynchTask.h.

Referenced by analyze(), and beginRun().