CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
27 
29 
31 
35 
36 #include "TMath.h"
37 #include <cmath>
38 
39 using namespace std;
40 using namespace edm;
41 
43  : // Get the debug parameter for verbose output
44  debug(pSet.getUntrackedParameter<bool>("debug", false)),
45  nMinHitsPhi(pSet.getUntrackedParameter<int>("nMinHitsPhi")),
46  maxAnglePhiSegm(pSet.getUntrackedParameter<double>("maxAnglePhiSegm")),
47  dt4DSegmentsToken_(consumes<DTRecSegment4DCollection>(pSet.getParameter<InputTag>("recoSegments"))),
48  muonGeomToken_(esConsumes<edm::Transition::BeginRun>()),
49  readLegacyVDriftDB(pSet.getParameter<bool>("readLegacyVDriftDB")) {
50  if (readLegacyVDriftDB) {
52  } else {
54  }
55 }
56 
58  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") << "DTRunConditionVar: destructor called";
59 
60  // free memory
61 }
62 
64  edm::Run const& iRun,
65  edm::EventSetup const& /* iSetup */) {
66  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar") << "DTRunConditionVar: bookHistograms";
67 
68  for (int wheel = -2; wheel <= 2; wheel++) {
69  for (int sec = 1; sec <= 14; sec++) {
70  for (int stat = 1; stat <= 4; stat++) {
71  bookChamberHistos(ibooker, DTChamberId(wheel, stat, sec), "VDrift_FromSegm", 100, 0.0043, 0.0065);
72  bookChamberHistos(ibooker, DTChamberId(wheel, stat, sec), "T0_FromSegm", 100, -25., 25.);
73  }
74  }
75  }
76 
77  return;
78 }
79 
81  // Get the DT Geometry
82  dtGeom = &setup.getData(muonGeomToken_);
83  return;
84 }
85 
86 void DTRunConditionVar::analyze(const Event& event, const EventSetup& eventSetup) {
87  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
88  << "--- [DTRunConditionVar] Event analysed #Run: " << event.id().run() << " #Event: " << event.id().event()
89  << endl;
90 
91  // Get the map of vdrift from the setup
92  if (readLegacyVDriftDB) {
93  mTimeMap_ = &eventSetup.getData(mTimeToken_);
94  vDriftMap_ = nullptr;
95  } else {
96  vDriftMap_ = &eventSetup.getData(vDriftToken_);
97  mTimeMap_ = nullptr;
98  // Consistency check: no parametrization is implemented for the time being
99  int version = vDriftMap_->version();
100  if (version != 1) {
101  throw cms::Exception("Configuration") << "only version 1 is presently supported for VDriftDB";
102  }
103  }
104 
105  // Get the segment collection from the event
106  Handle<DTRecSegment4DCollection> all4DSegments;
107  event.getByToken(dt4DSegmentsToken_, all4DSegments);
108 
109  // Loop over the segments
110  for (DTRecSegment4DCollection::const_iterator segment = all4DSegments->begin(); segment != all4DSegments->end();
111  ++segment) {
112  // Get the chamber from the setup
113  DTChamberId DTid = (DTChamberId)segment->chamberId();
114  uint32_t indexCh = DTid.rawId();
115 
116  // Fill v-drift values
117  if ((*segment).hasPhi()) {
118  int nHitsPhi = (*segment).phiSegment()->degreesOfFreedom() + 2;
119  double xdir = (*segment).phiSegment()->localDirection().x();
120  double zdir = (*segment).phiSegment()->localDirection().z();
121 
122  double anglePhiSegm = fabs(atan(xdir / zdir)) * 180. / TMath::Pi();
123 
124  if (nHitsPhi >= nMinHitsPhi && anglePhiSegm <= maxAnglePhiSegm) {
125  double segmentVDrift = segment->phiSegment()->vDrift();
126 
127  DTSuperLayerId indexSLPhi1(DTid, 1);
128  DTSuperLayerId indexSLPhi2(DTid, 3);
129 
130  float vDriftPhi1(0.), vDriftPhi2(0.);
131  float ResPhi1(0.), ResPhi2(0.);
132  if (readLegacyVDriftDB) { // Legacy format
133  int status1 = mTimeMap_->get(indexSLPhi1, vDriftPhi1, ResPhi1, DTVelocityUnits::cm_per_ns);
134  int status2 = mTimeMap_->get(indexSLPhi2, vDriftPhi2, ResPhi2, DTVelocityUnits::cm_per_ns);
135 
136  if (status1 != 0 || status2 != 0) {
137  DTSuperLayerId sl = (status1 != 0) ? indexSLPhi1 : indexSLPhi2;
138  throw cms::Exception("DTRunConditionVarClient") << "Could not find vDrift entry in DB for" << sl << endl;
139  }
140  } else {
141  vDriftPhi1 = vDriftMap_->get(DTWireId(indexSLPhi1.rawId()));
142  vDriftPhi2 = vDriftMap_->get(DTWireId(indexSLPhi2.rawId()));
143  }
144 
145  float vDriftMed = (vDriftPhi1 + vDriftPhi2) / 2.;
146 
147  segmentVDrift = vDriftMed * (1. - segmentVDrift);
148 
149  double segmentT0 = segment->phiSegment()->t0();
150 
151  if (segment->phiSegment()->ist0Valid())
152  (chamberHistos[indexCh])["T0_FromSegm"]->Fill(segmentT0);
153  if (segmentVDrift != vDriftMed)
154  (chamberHistos[indexCh])["VDrift_FromSegm"]->Fill(segmentVDrift);
155  }
156  }
157 
158  } //end loop on segment
159 
160 } //end analyze
161 
163  DQMStore::IBooker& ibooker, const DTChamberId& dtCh, string histoType, int nbins, float min, float max) {
164  int wh = dtCh.wheel();
165  int sc = dtCh.sector();
166  int st = dtCh.station();
167  stringstream wheel;
168  wheel << wh;
169  stringstream station;
170  station << st;
171  stringstream sector;
172  sector << sc;
173 
174  string bookingFolder = "DT/02-Segments/Wheel" + wheel.str() + "/Sector" + sector.str() + "/Station" + station.str();
175  string histoTag = "_W" + wheel.str() + "_Sec" + sector.str() + "_St" + station.str();
176 
177  ibooker.setCurrentFolder(bookingFolder);
178 
179  LogTrace("DTDQM|DTMonitorModule|DTRunConditionVar")
180  << "[DTRunConditionVar]: booking histos in " << bookingFolder << endl;
181 
182  string histoName = histoType + histoTag;
183  const string& histoLabel = histoType;
184 
185  (chamberHistos[dtCh.rawId()])[histoType] = ibooker.book1D(histoName, histoLabel, nbins, min, max);
186 }
187 
188 // Local Variables:
189 // show-trailing-whitespace: t
190 // truncate-lines: t
191 // End:
const double Pi
edm::EDGetTokenT< DTRecSegment4DCollection > dt4DSegmentsToken_
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
int version() const
Version numer specifying the structure of the payload. See .cc file for details.
constexpr uint32_t rawId() const
get the raw id
Definition: DetId.h:57
edm::ESGetToken< DTRecoConditions, DTRecoConditionsVdriftRcd > vDriftToken_
void analyze(const edm::Event &event, const edm::EventSetup &eventSetup) override
#define LogTrace(id)
const DTMtime * mTimeMap_
bool getData(T &iHolder) const
Definition: EventSetup.h:128
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:43
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
edm::ESGetToken< DTGeometry, MuonGeometryRecord > muonGeomToken_
void dqmBeginRun(const edm::Run &, const edm::EventSetup &) override
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
Transition
Definition: Transition.h:12
const DTRecoConditions * vDriftMap_
T min(T a, T b)
Definition: MathUtil.h:58
std::map< uint32_t, std::map< std::string, MonitorElement * > > chamberHistos
float get(const DTWireId &wireid, double *x=nullptr) const
Get the value correspoding to the given WireId, / using x[] as parameters of the parametrization whe...
#define debug
Definition: HDRShower.cc:19
int get(int wheelId, int stationId, int sectorId, int slId, float &mTime, float &mTrms, DTTimeUnits::type unit) const
Definition: DTMtime.cc:56
edm::ESGetToken< DTMtime, DTMtimeRcd > mTimeToken_
int sector() const
Definition: DTChamberId.h:49
~DTRunConditionVar() override
void bookChamberHistos(DQMStore::IBooker &, const DTChamberId &dtCh, std::string histoType, int, float, float)
const DTGeometry * dtGeom
MonitorElement * book1D(TString const &name, TString const &title, int const nchX, double const lowX, double const highX, FUNC onbooking=NOOP())
Definition: DQMStore.h:98
int station() const
Return the station number.
Definition: DTChamberId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
Definition: Run.h:45
DTRunConditionVar(const edm::ParameterSet &pset)