CMS 3D CMS Logo

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

#include <DTTTrigOffsetCalibration.h>

Inheritance diagram for DTTTrigOffsetCalibration:
edm::EDAnalyzer

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup)
 
void beginRun (const edm::Run &run, const edm::EventSetup &setup)
 
void bookHistos (DTChamberId)
 
 DTTTrigOffsetCalibration (const edm::ParameterSet &pset)
 
void endJob ()
 
virtual ~DTTTrigOffsetCalibration ()
 
- Public Member Functions inherited from edm::EDAnalyzer
 EDAnalyzer ()
 
std::string workerType () const
 
virtual ~EDAnalyzer ()
 

Private Types

typedef std::map< DTChamberId,
std::vector< TH1F * > > 
ChamberHistosMap
 

Private Attributes

bool checkNoisyChannels_
 
std::string dbLabel
 
bool doTTrigCorrection_
 
std::string theCalibChamber_
 
TFile * theFile_
 
double theMaxChi2_
 
double theMaxPhiAngle_
 
double theMaxZAngle_
 
edm::InputTag theRecHits4DLabel_
 
ChamberHistosMap theT0SegHistoMap_
 
const DTTtrigtTrigMap
 

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)
 
- Protected Member Functions inherited from edm::EDAnalyzer
CurrentProcessingContext const * currentContext () const
 

Detailed Description

No description available.

Date:
2010/02/16 10:03:23
Revision:
1.2
Author
A. Vilela Pereira

Definition at line 28 of file DTTTrigOffsetCalibration.h.

Member Typedef Documentation

typedef std::map<DTChamberId, std::vector<TH1F*> > DTTTrigOffsetCalibration::ChamberHistosMap
private

Definition at line 45 of file DTTTrigOffsetCalibration.h.

Constructor & Destructor Documentation

DTTTrigOffsetCalibration::DTTTrigOffsetCalibration ( const edm::ParameterSet pset)

Definition at line 41 of file DTTTrigOffsetCalibration.cc.

References edm::ParameterSet::getParameter(), and edm::ParameterSet::getUntrackedParameter().

41  {
42 
43  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Constructor called!";
44 
45  // the root file which will contain the histos
46  string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0SegHistos.root");
47  theFile_ = new TFile(rootFileName.c_str(), "RECREATE");
48  theFile_->cd();
49 
50  dbLabel = pset.getUntrackedParameter<string>("dbLabel", "");
51 
52  // Do t0-seg correction to ttrig
53  doTTrigCorrection_ = pset.getUntrackedParameter<bool>("doT0SegCorrection", false);
54 
55  // Chamber/s to calibrate
56  theCalibChamber_ = pset.getUntrackedParameter<string>("calibChamber", "All");
57 
58  // the name of the 4D rec hits collection
59  theRecHits4DLabel_ = pset.getParameter<InputTag>("recHits4DLabel");
60 
61  //get the switch to check the noisy channels
62  checkNoisyChannels_ = pset.getParameter<bool>("checkNoisyChannels");
63 
64  // get maximum chi2 value
65  theMaxChi2_ = pset.getParameter<double>("maxChi2");
66 
67  // Maximum incidence angle for Phi SL
68  theMaxPhiAngle_ = pset.getParameter<double>("maxAnglePhi");
69 
70  // Maximum incidence angle for Theta SL
71  theMaxZAngle_ = pset.getParameter<double>("maxAngleZ");
72 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTTTrigOffsetCalibration::~DTTTrigOffsetCalibration ( )
virtual

Definition at line 83 of file DTTTrigOffsetCalibration.cc.

83  {
84  theFile_->Close();
85  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Destructor called!";
86 }

Member Function Documentation

void DTTTrigOffsetCalibration::analyze ( const edm::Event event,
const edm::EventSetup eventSetup 
)
virtual

Implements edm::EDAnalyzer.

Definition at line 88 of file DTTTrigOffsetCalibration.cc.

References bookHistos(), DTChamberId, edm::EventSetup::get(), DTRecSegment2D::localDirection(), DTRecSegment2D::localPosition(), LogTrace, Geom::pi(), DTRecSegment2D::specificRecHits(), crabStatusFromReport::statusMap, DTChamber::superLayer(), DTSLRecSegment2D::superLayerId(), DTLayerId::superlayerId(), GeomDet::toGlobal(), GeomDet::toLocal(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().

88  {
89  theFile_->cd();
90  DTChamberId chosenChamberId;
91 
92  if(theCalibChamber_ != "All") {
93  stringstream linestr;
94  int selWheel, selStation, selSector;
95  linestr << theCalibChamber_;
96  linestr >> selWheel >> selStation >> selSector;
97  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
98  LogVerbatim("Calibration") << " chosen chamber " << chosenChamberId << endl;
99  }
100 
101  // Get the DT Geometry
102  ESHandle<DTGeometry> dtGeom;
103  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
104 
105  // Get the rechit collection from the event
106  Handle<DTRecSegment4DCollection> all4DSegments;
107  event.getByLabel(theRecHits4DLabel_, all4DSegments);
108 
109  // Get the map of noisy channels
111  if(checkNoisyChannels_) eventSetup.get<DTStatusFlagRcd>().get(statusMap);
112 
113  // Loop over segments by chamber
115  for (chamberIdIt = all4DSegments->id_begin();
116  chamberIdIt != all4DSegments->id_end();
117  ++chamberIdIt){
118 
119  // Get the chamber from the setup
120  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
121  LogTrace("Calibration") << "Chamber Id: " << *chamberIdIt;
122 
123  // Book histos
124  if(theT0SegHistoMap_.find(*chamberIdIt) == theT0SegHistoMap_.end()){
125  bookHistos(*chamberIdIt);
126  }
127 
128  // Calibrate just the chosen chamber/s
129  if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue;
130 
131  // Get the range for the corresponding ChamberId
132  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
133 
134  // Loop over the rechits of this DetUnit
135  for (DTRecSegment4DCollection::const_iterator segment = range.first;
136  segment!=range.second; ++segment){
137 
138  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
139  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
140 
141  //get the segment chi2
142  double chiSquare = ((*segment).chi2()/(*segment).degreesOfFreedom());
143  // cut on the segment chi2
144  if(chiSquare > theMaxChi2_) continue;
145 
146  // get the Phi 2D segment and plot the angle in the chamber RF
147  if(!((*segment).phiSegment())){
148  LogTrace("Calibration") << "No phi segment";
149  }
150  LocalPoint phiSeg2DPosInCham;
151  LocalVector phiSeg2DDirInCham;
152 
153  bool segmNoisy = false;
154  map<DTSuperLayerId,vector<DTRecHit1D> > hitsBySLMap;
155 
156  if((*segment).hasPhi()){
157  const DTChamberRecSegment2D* phiSeg = (*segment).phiSegment(); // phiSeg lives in the chamber RF
158  phiSeg2DPosInCham = phiSeg->localPosition();
159  phiSeg2DDirInCham = phiSeg->localDirection();
160 
161  vector<DTRecHit1D> phiHits = phiSeg->specificRecHits();
162  for(vector<DTRecHit1D>::const_iterator hit = phiHits.begin();
163  hit != phiHits.end(); ++hit) {
164  DTWireId wireId = (*hit).wireId();
165  DTSuperLayerId slId = wireId.superlayerId();
166  hitsBySLMap[slId].push_back(*hit);
167 
168  // Check for noisy channels to skip them
169  if(checkNoisyChannels_) {
170  bool isNoisy = false;
171  bool isFEMasked = false;
172  bool isTDCMasked = false;
173  bool isTrigMask = false;
174  bool isDead = false;
175  bool isNohv = false;
176  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
177  if(isNoisy) {
178  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
179  segmNoisy = true;
180  }
181  }
182  }
183  }
184 
185  // get the Theta 2D segment and plot the angle in the chamber RF
186  LocalVector zSeg2DDirInCham;
187  LocalPoint zSeg2DPosInCham;
188  if((*segment).hasZed()) {
189  const DTSLRecSegment2D* zSeg = (*segment).zSegment(); // zSeg lives in the SL RF
190  const DTSuperLayer* sl = chamber->superLayer(zSeg->superLayerId());
191  zSeg2DPosInCham = chamber->toLocal(sl->toGlobal((*zSeg).localPosition()));
192  zSeg2DDirInCham = chamber->toLocal(sl->toGlobal((*zSeg).localDirection()));
193  hitsBySLMap[zSeg->superLayerId()] = zSeg->specificRecHits();
194 
195  // Check for noisy channels to skip them
196  vector<DTRecHit1D> zHits = zSeg->specificRecHits();
197  for(vector<DTRecHit1D>::const_iterator hit = zHits.begin();
198  hit != zHits.end(); ++hit) {
199  DTWireId wireId = (*hit).wireId();
200  if(checkNoisyChannels_) {
201  bool isNoisy = false;
202  bool isFEMasked = false;
203  bool isTDCMasked = false;
204  bool isTrigMask = false;
205  bool isDead = false;
206  bool isNohv = false;
207  statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
208  if(isNoisy) {
209  LogTrace("Calibration") << "Wire: " << wireId << " is noisy, skipping!";
210  segmNoisy = true;
211  }
212  }
213  }
214  }
215 
216  if (segmNoisy) continue;
217 
218  LocalPoint segment4DLocalPos = (*segment).localPosition();
219  LocalVector segment4DLocalDir = (*segment).localDirection();
220  if(fabs(atan(segment4DLocalDir.y()/segment4DLocalDir.z())* 180./Geom::pi()) > theMaxZAngle_) continue; // cut on the angle
221  if(fabs(atan(segment4DLocalDir.x()/segment4DLocalDir.z())* 180./Geom::pi()) > theMaxPhiAngle_) continue; // cut on the angle
222  // Fill t0-seg values
223  if((*segment).hasPhi()) {
224  //if((segment->phiSegment()->t0()) != 0.00){
225  if(segment->phiSegment()->ist0Valid()){
226  (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
227  }
228  }
229  if((*segment).hasZed()){
230  //if((segment->zSegment()->t0()) != 0.00){
231  if(segment->zSegment()->ist0Valid()){
232  (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
233  }
234  }
235 
236  // Fill t0-seg values
237  // if((*segment).hasPhi()) (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
238  // if((*segment).hasZed()) (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
239  //if((*segment).hasZed() && (*segment).hasPhi()) {}
240 
241  /*//loop over the segments
242  for(map<DTSuperLayerId,vector<DTRecHit1D> >::const_iterator slIdAndHits = hitsBySLMap.begin(); slIdAndHits != hitsBySLMap.end(); ++slIdAndHits) {
243  if (slIdAndHits->second.size() < 3) continue;
244  DTSuperLayerId slId = slIdAndHits->first;
245 
246  }*/
247  }
248  }
249 }
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:52
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:49
T y() const
Definition: PV3DBase.h:57
LocalPoint toLocal(const GlobalPoint &gp) const
Conversion to the R.F. of the GeomDet.
Definition: GeomDet.h:64
identifier iterator
Definition: RangeMap.h:138
DTSuperLayerId superlayerId() const
Return the corresponding SuperLayerId.
Definition: DTLayerId.h:61
C::const_iterator const_iterator
constant access iterator type
Definition: RangeMap.h:45
T z() const
Definition: PV3DBase.h:58
#define LogTrace(id)
std::vector< DTRecHit1D > specificRecHits() const
Access to specific components.
DTSuperLayerId superLayerId() const
The id of the superlayer on which reside the segment.
virtual LocalPoint localPosition() const
local position in SL frame
const T & get() const
Definition: EventSetup.h:55
double pi()
Definition: Pi.h:31
virtual LocalVector localDirection() const
the local direction in SL frame
T x() const
Definition: PV3DBase.h:56
const DTSuperLayer * superLayer(DTSuperLayerId id) const
Return the superlayer corresponding to the given id.
Definition: DTChamber.cc:67
void DTTTrigOffsetCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 74 of file DTTTrigOffsetCalibration.cc.

References edm::EventSetup::get().

74  {
76  ESHandle<DTTtrig> tTrig;
77  setup.get<DTTtrigRcd>().get(dbLabel,tTrig);
78  tTrigMap = &*tTrig;
79  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]: TTrig version: " << tTrig->version() << endl;
80  }
81 }
const T & get() const
Definition: EventSetup.h:55
void DTTTrigOffsetCalibration::bookHistos ( DTChamberId  chId)

Definition at line 298 of file DTTTrigOffsetCalibration.cc.

References combine::histos, LogTrace, DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, and DTChamberId::wheel().

298  {
299 
300  LogTrace("Calibration") << " Booking histos for Chamber: " << chId;
301 
302  // Compose the chamber name
303  stringstream wheel; wheel << chId.wheel();
304  stringstream station; station << chId.station();
305  stringstream sector; sector << chId.sector();
306 
307  string chHistoName =
308  "_W" + wheel.str() +
309  "_St" + station.str() +
310  "_Sec" + sector.str();
311 
312  /*// Define the step
313  stringstream Step; Step << step;
314 
315  string chHistoName =
316  "_STEP" + Step.str() +
317  "_W" + wheel.str() +
318  "_St" + station.str() +
319  "_Sec" + sector.str();
320 
321  theDbe->setCurrentFolder("DT/DTCalibValidation/Wheel" + wheel.str() +
322  "/Station" + station.str() +
323  "/Sector" + sector.str());
324  // Create the monitor elements
325  vector<MonitorElement *> histos;
326  // Note hte order matters
327  histos.push_back(theDbe->book1D("hRPhiSegT0"+chHistoName, "t0 from Phi segments", 200, -25., 25.));
328  histos.push_back(theDbe->book1D("hRZSegT0"+chHistoName, "t0 from Z segments", 200, -25., 25.));*/
329 
330  vector<TH1F*> histos;
331  // Note the order matters
332  histos.push_back(new TH1F(("hRPhiSegT0"+chHistoName).c_str(), "t0 from Phi segments", 250, -60., 60.));
333  if(chId.station() != 4) histos.push_back(new TH1F(("hRZSegT0"+chHistoName).c_str(), "t0 from Z segments", 250, -60., 60.));
334 
335  theT0SegHistoMap_[chId] = histos;
336 }
#define LogTrace(id)
dictionary histos
Definition: combine.py:3
int sector() const
Definition: DTChamberId.h:63
int station() const
Return the station number.
Definition: DTChamberId.h:53
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:47
void DTTTrigOffsetCalibration::endJob ( void  )
virtual

Reimplemented from edm::EDAnalyzer.

Definition at line 251 of file DTTTrigOffsetCalibration.cc.

References DTSuperLayerId, DTTimeUnits::ns, DTTtrig::set(), DTChamberId::station(), and DTCalibDBUtils::writeToDB().

251  {
252  theFile_->cd();
253 
254  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]Writing histos to file!" << endl;
255 
256  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
257  for(vector<TH1F*>::const_iterator itHist = (*itChHistos).second.begin();
258  itHist != (*itChHistos).second.end(); ++itHist) (*itHist)->Write();
259  }
260 
261  if(doTTrigCorrection_){
262  // Create the object to be written to DB
263  DTTtrig* tTrig = new DTTtrig();
264 
265  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
266  DTChamberId chId = itChHistos->first;
267  // Get SuperLayerId's for each ChamberId
268  vector<DTSuperLayerId> slIds;
269  slIds.push_back(DTSuperLayerId(chId,1));
270  slIds.push_back(DTSuperLayerId(chId,3));
271  if(chId.station() != 4) slIds.push_back(DTSuperLayerId(chId,2));
272 
273  for(vector<DTSuperLayerId>::const_iterator itSl = slIds.begin(); itSl != slIds.end(); ++itSl){
274  // Get old values from DB
275  float ttrigMean = 0;
276  float ttrigSigma = 0;
277  float kFactor = 0;
278  tTrigMap->get(*itSl,ttrigMean,ttrigSigma,kFactor,DTTimeUnits::ns);
279  //FIXME: verify if values make sense
280  // Set new values
281  float ttrigMeanNew = ttrigMean;
282  float ttrigSigmaNew = ttrigSigma;
283  float t0SegMean = (itSl->superLayer() != 2)?itChHistos->second[0]->GetMean():itChHistos->second[1]->GetMean();
284 
285  float kFactorNew = (kFactor*ttrigSigma+t0SegMean)/ttrigSigma;
286 
287  tTrig->set(*itSl,ttrigMeanNew,ttrigSigmaNew,kFactorNew,DTTimeUnits::ns);
288  }
289  }
290  LogVerbatim("Calibration")<< "[DTTTrigOffsetCalibration]Writing ttrig object to DB!" << endl;
291  // Write the object to DB
292  string tTrigRecord = "DTTtrigRcd";
293  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
294  }
295 }
int set(int wheelId, int stationId, int sectorId, int slId, float tTrig, float tTrms, float kFact, DTTimeUnits::type unit)
Definition: DTTtrig.cc:262
DTSuperLayerId
int get(int wheelId, int stationId, int sectorId, int slId, float &tTrig, float &tTrms, float &kFact, DTTimeUnits::type unit) const
get content
Definition: DTTtrig.cc:87
int station() const
Return the station number.
Definition: DTChamberId.h:53
static void writeToDB(std::string record, T *payload)

Member Data Documentation

bool DTTTrigOffsetCalibration::checkNoisyChannels_
private

Definition at line 60 of file DTTTrigOffsetCalibration.h.

std::string DTTTrigOffsetCalibration::dbLabel
private

Definition at line 77 of file DTTTrigOffsetCalibration.h.

bool DTTTrigOffsetCalibration::doTTrigCorrection_
private

Definition at line 54 of file DTTTrigOffsetCalibration.h.

std::string DTTTrigOffsetCalibration::theCalibChamber_
private

Definition at line 75 of file DTTTrigOffsetCalibration.h.

TFile* DTTTrigOffsetCalibration::theFile_
private

Definition at line 51 of file DTTTrigOffsetCalibration.h.

double DTTTrigOffsetCalibration::theMaxChi2_
private

Definition at line 66 of file DTTTrigOffsetCalibration.h.

double DTTTrigOffsetCalibration::theMaxPhiAngle_
private

Definition at line 69 of file DTTTrigOffsetCalibration.h.

double DTTTrigOffsetCalibration::theMaxZAngle_
private

Definition at line 72 of file DTTTrigOffsetCalibration.h.

edm::InputTag DTTTrigOffsetCalibration::theRecHits4DLabel_
private

Definition at line 48 of file DTTTrigOffsetCalibration.h.

ChamberHistosMap DTTTrigOffsetCalibration::theT0SegHistoMap_
private

Definition at line 63 of file DTTTrigOffsetCalibration.h.

const DTTtrig* DTTTrigOffsetCalibration::tTrigMap
private

Definition at line 57 of file DTTTrigOffsetCalibration.h.