CMS 3D CMS Logo

DTRunConditionVar.cc
Go to the documentation of this file.
1 
2 /******* \class DTRunConditionVar *******
3  *
4  * Description:
5  *
6  * detailed description
7  *
8  * \author : Paolo Bellan, Antonio Branca
9  * $date : 23/09/2011 15:42:04 CET $
10  *
11  * Modification:
12  *
13  *********************************/
14 
15 #include "DTRunConditionVar.h"
16 
21 
23 
25 
28 
31 
33 
35 
39 
40 #include <TMath.h>
41 #include <cmath>
42 
43 using namespace std;
44 using namespace edm;
45 
47  // Get the debug parameter for verbose output
48  debug(pSet.getUntrackedParameter<bool>("debug",false)),
49  nMinHitsPhi(pSet.getUntrackedParameter<int>("nMinHitsPhi")),
50  maxAnglePhiSegm(pSet.getUntrackedParameter<double>("maxAnglePhiSegm")),
51  dt4DSegmentsToken_(consumes<DTRecSegment4DCollection>(
52  pSet.getParameter<InputTag>("recoSegments")))
53 {
54 
55 }
56 
58 {
59  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
60  << "DTRunConditionVar: destructor called";
61 
62  // free memory
63 }
64 
66  edm::Run const & iRun,
67  edm::EventSetup const & /* iSetup */) {
68 
69  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
70  << "DTRunConditionVar: bookHistograms";
71 
72  for(int wheel=-2;wheel<=2;wheel++){
73  for(int sec=1; sec<=14; sec++) {
74  for(int stat=1; stat<=4; stat++) {
75 
76  bookChamberHistos(ibooker,DTChamberId(wheel,stat,sec),"VDrift_FromSegm",100,0.0043,0.0065);
77  bookChamberHistos(ibooker,DTChamberId(wheel,stat,sec),"T0_FromSegm",100,-25.,25.);
78 
79  }
80  }
81  }
82 
83 
84  return;
85 }
86 
88 {
89  // Get the DT Geometry
90  setup.get<MuonGeometryRecord>().get(dtGeom);
91 
92  return;
93 }
94 
96  const EventSetup& eventSetup)
97 {
98 
99  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") <<
100  "--- [DTRunConditionVar] Event analysed #Run: " <<
101  event.id().run() << " #Event: " << event.id().event() << endl;
102 
103  // Get the DT Geometry
105  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
106 
107  // Get the map of vdrift from the setup
108  eventSetup.get<DTMtimeRcd>().get(mTime);
109  mTimeMap_ = &*mTime;
110 
111  // Get the segment collection from the event
112  Handle<DTRecSegment4DCollection> all4DSegments;
113  event.getByToken(dt4DSegmentsToken_, all4DSegments);
114 
115  // Loop over the segments
116  for(DTRecSegment4DCollection::const_iterator segment = all4DSegments->begin();
117  segment != all4DSegments->end(); ++segment){
118 
119  // Get the chamber from the setup
120  DTChamberId DTid = (DTChamberId) segment->chamberId();
121  uint32_t indexCh = DTid.rawId();
122 
123  // Fill v-drift values
124  if( (*segment).hasPhi() ) {
125 
126  int nHitsPhi = (*segment).phiSegment()->degreesOfFreedom()+2;
127  double xdir = (*segment).phiSegment()->localDirection().x();
128  double zdir = (*segment).phiSegment()->localDirection().z();
129 
130  double anglePhiSegm = fabs(atan(xdir/zdir))*180./TMath::Pi();
131 
132  if( nHitsPhi >= nMinHitsPhi && anglePhiSegm <= maxAnglePhiSegm ) {
133 
134  double segmentVDrift = segment->phiSegment()->vDrift();
135 
136 
137  DTSuperLayerId indexSLPhi1(DTid,1);
138  DTSuperLayerId indexSLPhi2(DTid,3);
139 
140  float vDriftPhi1(0.), vDriftPhi2(0.);
141  float ResPhi1(0.), ResPhi2(0.);
142  int status1 = mTimeMap_->get(indexSLPhi1,vDriftPhi1,ResPhi1,DTVelocityUnits::cm_per_ns);
143  int status2 = mTimeMap_->get(indexSLPhi2,vDriftPhi2,ResPhi2,DTVelocityUnits::cm_per_ns);
144 
145  if(status1 != 0 || status2 != 0) {
146  DTSuperLayerId sl = (status1 != 0) ? indexSLPhi1 : indexSLPhi2;
147  throw cms::Exception("DTRunConditionVarClient") << "Could not find vDrift entry in DB for"
148  << sl << endl;
149  }
150 
151  float vDriftMed = (vDriftPhi1 + vDriftPhi2) / 2.;
152 
153  segmentVDrift = vDriftMed*(1. - segmentVDrift);
154 
155  double segmentT0 = segment->phiSegment()->t0();
156 
157  if( segment->phiSegment()->ist0Valid() ) (chamberHistos[indexCh])["T0_FromSegm"]->Fill(segmentT0);
158  if( segmentVDrift != vDriftMed ) (chamberHistos[indexCh])["VDrift_FromSegm"]->Fill(segmentVDrift);
159 
160  }
161  }
162 
163  } //end loop on segment
164 
165 } //end analyze
166 
167 void DTRunConditionVar::bookChamberHistos(DQMStore::IBooker & ibooker,const DTChamberId& dtCh, string histoType, int nbins, float min, float max) {
168 
169  int wh = dtCh.wheel();
170  int sc = dtCh.sector();
171  int st = dtCh.station();
172  stringstream wheel; wheel << wh;
173  stringstream station; station << st;
174  stringstream sector; sector << sc;
175 
176  string bookingFolder = "DT/02-Segments/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str();
177  string histoTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
178 
179  ibooker.setCurrentFolder(bookingFolder);
180 
181  LogTrace ("DTDQM|DTMonitorModule|DTRunConditionVar")
182  << "[DTRunConditionVar]: booking histos in " << bookingFolder << endl;
183 
184  string histoName = histoType + histoTag;
185  string histoLabel = histoType;
186 
187  (chamberHistos[dtCh.rawId()])[histoType] =
188  ibooker.book1D(histoName,histoLabel,nbins,min,max);
189 
190 }
191 
192 // Local Variables:
193 // show-trailing-whitespace: t
194 // truncate-lines: t
195 // End:
const double Pi
edm::ESHandle< DTMtime > mTime
edm::EDGetTokenT< DTRecSegment4DCollection > dt4DSegmentsToken_
def setup(process, global_tag, zero_tesla=False)
Definition: GeneralSetup.py:1
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
const DTMtime * mTimeMap_
uint32_t rawId() const
get the raw id
Definition: DetId.h:44
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:118
T min(T a, T b)
Definition: MathUtil.h:58
#define LogTrace(id)
std::map< uint32_t, std::map< std::string, MonitorElement * > > chamberHistos
#define debug
Definition: HDRShower.cc:19
void setCurrentFolder(const std::string &fullpath)
Definition: DQMStore.cc:279
const T & get() const
Definition: EventSetup.h:59
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:82
HLT enums.
int sector() const
Definition: DTChamberId.h:61
~DTRunConditionVar() override
void bookChamberHistos(DQMStore::IBooker &, const DTChamberId &dtCh, std::string histoType, int, float, float)
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
edm::ESHandle< DTGeometry > dtGeom
Definition: event.py:1
Definition: Run.h:43
DTRunConditionVar(const edm::ParameterSet &pset)