CMS 3D CMS Logo

DTTPAnalyzer.cc
Go to the documentation of this file.
1 
12 
13 #include <string>
14 #include <map>
15 
16 //class DTT0;
17 class DTLayerId;
18 class DTWireId;
19 class DTGeometry;
20 class DTTTrigBaseSync;
21 class TFile;
22 
24 public:
26  ~DTTPAnalyzer() override;
27 
28  void analyze(const edm::Event&, const edm::EventSetup&) override;
29  void endJob() override;
30 
31 private:
33 
34  const bool subtractT0_;
36 
37  TFile* rootFile_;
40  std::unique_ptr<DTTTrigBaseSync> tTrigSync_;
41 
42  // Map of the t0 and sigma histos by layer
43  std::map<DTWireId, int> nDigisPerWire_;
44  std::map<DTWireId, double> sumWPerWire_;
45  std::map<DTWireId, double> sumW2PerWire_;
46 };
47 
51 
54 
57 
59 
61 
62 #include "TH1F.h"
63 #include "TFile.h"
64 
66  : subtractT0_(pset.getParameter<bool>("subtractT0")),
67  digiToken_(consumes<DTDigiCollection>(pset.getParameter<edm::InputTag>("digiLabel"))),
68  dtGeomToken_(esConsumes()) {
69  std::string rootFileName = pset.getUntrackedParameter<std::string>("rootFileName");
70  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
71  rootFile_->cd();
72 
73  if (subtractT0_)
74  tTrigSync_ = DTTTrigSyncFactory::get()->create(pset.getParameter<std::string>("tTrigMode"),
75  pset.getParameter<edm::ParameterSet>("tTrigModeConfig"),
77 }
78 
80 
82  // Get the digis from the event
83  const edm::Handle<DTDigiCollection>& digis = event.getHandle(digiToken_);
84 
85  if (subtractT0_) {
86  tTrigSync_->setES(setup);
87  }
88  // Get the DT Geometry
89  dtGeom_ = setup.getHandle(dtGeomToken_);
90 
91  // Iterate through all digi collections ordered by LayerId
93  for (dtLayerIt = digis->begin(); dtLayerIt != digis->end(); ++dtLayerIt) {
94  // Get the iterators over the digis associated with this LayerId
95  const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;
96 
97  // Get the layerId
98  const DTLayerId layerId = (*dtLayerIt).first; //FIXME: check to be in the right sector
99 
100  // Loop over all digis in the given layer
101  for (DTDigiCollection::const_iterator digi = digiRange.first; digi != digiRange.second; ++digi) {
102  const DTWireId wireId(layerId, (*digi).wire());
103 
104  double t0 = (*digi).countsTDC();
105 
106  //FIXME: Reject digis not coming from TP
107 
108  if (subtractT0_) {
109  const DTLayer* layer = nullptr; //fake
110  const GlobalPoint glPt; //fake
111  double offset = tTrigSync_->offset(layer, wireId, glPt);
112  t0 -= offset;
113  }
114 
115  if (nDigisPerWire_.find(wireId) == nDigisPerWire_.end()) {
116  nDigisPerWire_[wireId] = 0;
117  sumWPerWire_[wireId] = 0.;
118  sumW2PerWire_[wireId] = 0.;
119  }
120 
121  ++nDigisPerWire_[wireId];
122  sumWPerWire_[wireId] += t0;
123  sumW2PerWire_[wireId] += t0 * t0;
124  }
125  }
126 }
127 
129  rootFile_->cd();
130  std::map<DTLayerId, TH1F*> meanHistoMap;
131  std::map<DTLayerId, TH1F*> sigmaHistoMap;
132  for (std::map<DTWireId, int>::const_iterator wireIdIt = nDigisPerWire_.begin(); wireIdIt != nDigisPerWire_.end();
133  ++wireIdIt) {
134  DTWireId wireId((*wireIdIt).first);
135 
136  int nDigis = nDigisPerWire_[wireId];
137  double sumW = sumWPerWire_[wireId];
138  double sumW2 = sumW2PerWire_[wireId];
139 
140  double mean = sumW / nDigis;
141  double rms = sumW2 / nDigis - mean * mean;
142  rms = sqrt(rms);
143 
144  DTLayerId layerId = wireId.layerId();
145  if (meanHistoMap.find(layerId) == meanHistoMap.end()) {
147  const int firstChannel = dtGeom_->layer(layerId)->specificTopology().firstChannel();
148  const int nWires = dtGeom_->layer(layerId)->specificTopology().channels();
149  TH1F* meanHistoTP = new TH1F((histoName + "_tpMean").c_str(),
150  "mean from test pulses by channel",
151  nWires,
152  firstChannel,
153  (firstChannel + nWires));
154  TH1F* sigmaHistoTP = new TH1F((histoName + "_tpSigma").c_str(),
155  "sigma from test pulses by channel",
156  nWires,
157  firstChannel,
158  (firstChannel + nWires));
159  meanHistoMap[layerId] = meanHistoTP;
160  sigmaHistoMap[layerId] = sigmaHistoTP;
161  }
162  // Fill the histograms
163  int nBin = meanHistoMap[layerId]->GetXaxis()->FindFixBin(wireId.wire());
164  meanHistoMap[layerId]->SetBinContent(nBin, mean);
165  sigmaHistoMap[layerId]->SetBinContent(nBin, rms);
166  }
167 
168  for (std::map<DTLayerId, TH1F*>::const_iterator key = meanHistoMap.begin(); key != meanHistoMap.end(); ++key) {
169  meanHistoMap[(*key).first]->Write();
170  sigmaHistoMap[(*key).first]->Write();
171  }
172 }
173 
176  std::stringstream theStream;
177  theStream << "Ch_" << lId.wheel() << "_" << lId.station() << "_" << lId.sector() << "_SL" << lId.superlayer() << "_L"
178  << lId.layer();
179  theStream >> histoName;
180  return histoName;
181 }
182 
int station() const
Return the station number.
Definition: DTChamberId.h:42
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
int wire() const
Return the wire number.
Definition: DTWireId.h:42
const edm::ESGetToken< DTGeometry, MuonGeometryRecord > dtGeomToken_
Definition: DTTPAnalyzer.cc:39
DTTPAnalyzer(const edm::ParameterSet &)
Definition: DTTPAnalyzer.cc:65
std::map< DTWireId, double > sumWPerWire_
Definition: DTTPAnalyzer.cc:44
const bool subtractT0_
Definition: DTTPAnalyzer.cc:34
std::map< DTWireId, int > nDigisPerWire_
Definition: DTTPAnalyzer.cc:43
TFile * rootFile_
Definition: DTTPAnalyzer.cc:37
std::string getHistoName(const DTLayerId &)
int firstChannel() const
Returns the wire number of the first wire.
Definition: DTTopology.h:79
void endJob() override
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
T sqrt(T t)
Definition: SSEVec.h:19
~DTTPAnalyzer() override
Definition: DTTPAnalyzer.cc:79
std::unique_ptr< DTTTrigBaseSync > tTrigSync_
Definition: DTTPAnalyzer.cc:40
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const DTTopology & specificTopology() const
Definition: DTLayer.cc:37
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: DTTPAnalyzer.cc:81
std::map< DTWireId, double > sumW2PerWire_
Definition: DTTPAnalyzer.cc:45
int superlayer() const
Return the superlayer number (deprecated method name)
std::pair< const_iterator, const_iterator > Range
std::vector< DigiType >::const_iterator const_iterator
int layer() const
Return the layer number.
Definition: DTLayerId.h:42
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:39
HLT enums.
int sector() const
Definition: DTChamberId.h:49
const edm::EDGetTokenT< DTDigiCollection > digiToken_
Definition: DTTPAnalyzer.cc:35
edm::ESHandle< DTGeometry > dtGeom_
Definition: DTTPAnalyzer.cc:38
#define get
DTLayerId layerId() const
Return the corresponding LayerId.
Definition: DTWireId.h:45
int channels() const
Returns the number of wires in the layer.
Definition: DTTopology.h:76
Definition: event.py:1
const DTLayer * layer(const DTLayerId &id) const
Return a layer given its id.
Definition: DTGeometry.cc:96