CMS 3D CMS Logo

List of all members | Public Member Functions | Private Types | Private Member Functions | Private Attributes
DTTTrigOffsetCalibration Class Reference

#include <DTTTrigOffsetCalibration.h>

Inheritance diagram for DTTTrigOffsetCalibration:
edm::EDAnalyzer edm::EDConsumerBase

Public Member Functions

void analyze (const edm::Event &event, const edm::EventSetup &eventSetup) override
 
void beginRun (const edm::Run &run, const edm::EventSetup &setup) override
 
 DTTTrigOffsetCalibration (const edm::ParameterSet &pset)
 
void endJob () override
 
 ~DTTTrigOffsetCalibration () override
 
- Public Member Functions inherited from edm::EDAnalyzer
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
 EDAnalyzer ()
 
SerialTaskQueueglobalLuminosityBlocksQueue ()
 
SerialTaskQueueglobalRunsQueue ()
 
ModuleDescription const & moduleDescription () const
 
std::string workerType () const
 
 ~EDAnalyzer () override
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
ProductResolverIndexAndSkipBit uncheckedIndexFrom (EDGetToken) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Types

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

Private Member Functions

void bookHistos (DTChamberId)
 

Private Attributes

std::string dbLabel_
 
bool doTTrigCorrection_
 
TFile * rootFile_
 
DTSegmentSelector select_
 
std::string theCalibChamber_
 
edm::InputTag theRecHits4DLabel_
 
ChamberHistosMap theT0SegHistoMap_
 
const DTTtrigtTrigMap_
 

Additional Inherited Members

- Public Types inherited from edm::EDAnalyzer
typedef EDAnalyzer ModuleType
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::EDAnalyzer
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
static bool wantsGlobalLuminosityBlocks ()
 
static bool wantsGlobalRuns ()
 
static bool wantsStreamLuminosityBlocks ()
 
static bool wantsStreamRuns ()
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes ()
 
template<typename ESProduct , typename ESRecord , Transition Tr = Transition::Event>
auto esConsumes (ESInputTag const &tag)
 
template<typename ESProduct , Transition Tr = Transition::Event>
auto esConsumes (eventsetup::EventSetupRecordKey const &, ESInputTag const &tag)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

No description available.

Author
A. Vilela Pereira

Definition at line 27 of file DTTTrigOffsetCalibration.h.

Member Typedef Documentation

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

Definition at line 39 of file DTTTrigOffsetCalibration.h.

Constructor & Destructor Documentation

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

Definition at line 37 of file DTTTrigOffsetCalibration.cc.

References edm::ParameterSet::getUntrackedParameter(), rootFile_, and DTAnalyzerDetailed_cfi::rootFileName.

37  :
38  select_(pset),
39  theRecHits4DLabel_(pset.getParameter<InputTag>("recHits4DLabel")),
40  doTTrigCorrection_(pset.getUntrackedParameter<bool>("doT0SegCorrection", false)),
41  theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")),
42  dbLabel_(pset.getUntrackedParameter<string>("dbLabel", "")) {
43 
44  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Constructor called!";
45 
46  // the root file which will contain the histos
47  string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0SegHistos.root");
48  rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
49  rootFile_->cd();
50 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
DTTTrigOffsetCalibration::~DTTTrigOffsetCalibration ( )
override

Definition at line 61 of file DTTTrigOffsetCalibration.cc.

References rootFile_.

61  {
62  rootFile_->Close();
63  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Destructor called!";
64 }

Member Function Documentation

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

Definition at line 66 of file DTTTrigOffsetCalibration.cc.

References bookHistos(), relativeConstraints::chamber, DTGeometry::chamber(), HcalObjRepresent::Fill(), edm::EventSetup::get(), LogTrace, rootFile_, select_, theCalibChamber_, theRecHits4DLabel_, theT0SegHistoMap_, and GeomDet::toGlobal().

66  {
67  rootFile_->cd();
68  DTChamberId chosenChamberId;
69 
70  if(theCalibChamber_ != "All") {
71  stringstream linestr;
72  int selWheel, selStation, selSector;
73  linestr << theCalibChamber_;
74  linestr >> selWheel >> selStation >> selSector;
75  chosenChamberId = DTChamberId(selWheel, selStation, selSector);
76  LogVerbatim("Calibration") << " chosen chamber " << chosenChamberId << endl;
77  }
78 
79  // Get the DT Geometry
80  ESHandle<DTGeometry> dtGeom;
81  eventSetup.get<MuonGeometryRecord>().get(dtGeom);
82 
83  // Get the rechit collection from the event
85  event.getByLabel(theRecHits4DLabel_, all4DSegments);
86 
87  // Loop over segments by chamber
88  DTRecSegment4DCollection::id_iterator chamberIdIt;
89  for(chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt){
90 
91  // Get the chamber from the setup
92  const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
93  LogTrace("Calibration") << "Chamber Id: " << *chamberIdIt;
94 
95  // Book histos
96  if(theT0SegHistoMap_.find(*chamberIdIt) == theT0SegHistoMap_.end()){
97  bookHistos(*chamberIdIt);
98  }
99 
100  // Calibrate just the chosen chamber/s
101  if((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId)) continue;
102 
103  // Get the range for the corresponding ChamberId
104  DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
105 
106  // Loop over the rechits of this DetUnit
107  for(DTRecSegment4DCollection::const_iterator segment = range.first;
108  segment != range.second; ++segment){
109 
110  LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
111  << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
112 
113  if( !select_(*segment, event, eventSetup) ) continue;
114 
115  // Fill t0-seg values
116  if( (*segment).hasPhi() ) {
117  //if( segment->phiSegment()->ist0Valid() ){
118  if( (segment->phiSegment()->t0()) != 0.00 ){
119  (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
120  }
121  }
122  if( (*segment).hasZed() ){
123  //if( segment->zSegment()->ist0Valid() ){
124  if( (segment->zSegment()->t0()) != 0.00 ){
125  (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
126  }
127  }
128  } // DTRecSegment4DCollection::const_iterator segment
129  } // DTRecSegment4DCollection::id_iterator chamberIdIt
130 } // DTTTrigOffsetCalibration::analyze
std::pair< const_iterator, const_iterator > range
iterator range
Definition: RangeMap.h:50
const DTChamber * chamber(const DTChamberId &id) const
Return a DTChamber given its id.
Definition: DTGeometry.cc:117
GlobalPoint toGlobal(const Local2DPoint &lp) const
Conversion to the global R.F. from the R.F. of the GeomDet.
Definition: GeomDet.h:54
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
#define LogTrace(id)
T get() const
Definition: EventSetup.h:68
void DTTTrigOffsetCalibration::beginRun ( const edm::Run run,
const edm::EventSetup setup 
)
override

Definition at line 52 of file DTTTrigOffsetCalibration.cc.

References dbLabel_, doTTrigCorrection_, edm::EventSetup::get(), tTrigMap_, and DTTtrig::version().

52  {
54  ESHandle<DTTtrig> tTrig;
55  setup.get<DTTtrigRcd>().get(dbLabel_,tTrig);
56  tTrigMap_ = &*tTrig;
57  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]: TTrig version: " << tTrig->version() << endl;
58  }
59 }
const std::string & version() const
access version
Definition: DTTtrig.cc:231
T get() const
Definition: EventSetup.h:68
void DTTTrigOffsetCalibration::bookHistos ( DTChamberId  chId)
private

Definition at line 179 of file DTTTrigOffsetCalibration.cc.

References plotFactory::histos, LogTrace, DTChamberId::sector(), DTChamberId::station(), relativeConstraints::station, AlCaHLTBitMon_QueryRunRegistry::string, theT0SegHistoMap_, DTChamberId::wheel(), and makeMuonMisalignmentScenario::wheel.

Referenced by analyze().

179  {
180 
181  LogTrace("Calibration") << " Booking histos for Chamber: " << chId;
182 
183  // Compose the chamber name
184  std::string wheel = std::to_string(chId.wheel());
185  std::string station = std::to_string(chId.station());
186  std::string sector = std::to_string(chId.sector());
187 
188  string chHistoName =
189  "_W" + wheel +
190  "_St" + station +
191  "_Sec" + sector;
192 
193  vector<TH1F*> histos;
194  // Note the order matters
195  histos.push_back(new TH1F(("hRPhiSegT0"+chHistoName).c_str(), "t0 from Phi segments", 500, -60., 60.));
196  if(chId.station() != 4) histos.push_back(new TH1F(("hRZSegT0"+chHistoName).c_str(), "t0 from Z segments", 500, -60., 60.));
197 
198  theT0SegHistoMap_[chId] = histos;
199 }
#define LogTrace(id)
int sector() const
Definition: DTChamberId.h:61
int station() const
Return the station number.
Definition: DTChamberId.h:51
int wheel() const
Return the wheel number.
Definition: DTChamberId.h:45
void DTTTrigOffsetCalibration::endJob ( void  )
overridevirtual

Reimplemented from edm::EDAnalyzer.

Definition at line 132 of file DTTTrigOffsetCalibration.cc.

References doTTrigCorrection_, DTTtrig::get(), DTAnalyzerDetailed_cfi::kFactor, DTTimeUnits::ns, rootFile_, DTTtrig::set(), DTChamberId::station(), theT0SegHistoMap_, tTrigMap_, and DTCalibDBUtils::writeToDB().

Referenced by o2olib.O2ORunMgr::executeJob().

132  {
133  rootFile_->cd();
134 
135  LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Writing histos to file!" << endl;
136 
137  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
138  for(vector<TH1F*>::const_iterator itHist = (*itChHistos).second.begin();
139  itHist != (*itChHistos).second.end(); ++itHist) (*itHist)->Write();
140  }
141 
142  if(doTTrigCorrection_){
143  // Create the object to be written to DB
144  DTTtrig* tTrig = new DTTtrig();
145 
146  for(ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end(); ++itChHistos){
147  DTChamberId chId = itChHistos->first;
148  // Get SuperLayerId's for each ChamberId
149  vector<DTSuperLayerId> slIds;
150  slIds.push_back(DTSuperLayerId(chId,1));
151  slIds.push_back(DTSuperLayerId(chId,3));
152  if(chId.station() != 4) slIds.push_back(DTSuperLayerId(chId,2));
153 
154  for(vector<DTSuperLayerId>::const_iterator itSl = slIds.begin(); itSl != slIds.end(); ++itSl){
155  // Get old values from DB
156  float ttrigMean = 0;
157  float ttrigSigma = 0;
158  float kFactor = 0;
159  tTrigMap_->get(*itSl,ttrigMean,ttrigSigma,kFactor,DTTimeUnits::ns);
160  //FIXME: verify if values make sense
161  // Set new values
162  float ttrigMeanNew = ttrigMean;
163  float ttrigSigmaNew = ttrigSigma;
164  float t0SegMean = (itSl->superLayer() != 2)?itChHistos->second[0]->GetMean():itChHistos->second[1]->GetMean();
165 
166  float kFactorNew = (kFactor*ttrigSigma+t0SegMean)/ttrigSigma;
167 
168  tTrig->set(*itSl,ttrigMeanNew,ttrigSigmaNew,kFactorNew,DTTimeUnits::ns);
169  }
170  }
171  LogVerbatim("Calibration")<< "[DTTTrigOffsetCalibration] Writing ttrig object to DB!" << endl;
172  // Write the object to DB
173  string tTrigRecord = "DTTtrigRcd";
174  DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
175  }
176 }
int set(int wheelId, int stationId, int sectorId, int slId, float tTrig, float tTrms, float kFact, DTTimeUnits::type unit)
Definition: DTTtrig.cc:248
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:85
int station() const
Return the station number.
Definition: DTChamberId.h:51
static void writeToDB(std::string record, T *payload)

Member Data Documentation

std::string DTTTrigOffsetCalibration::dbLabel_
private

Definition at line 47 of file DTTTrigOffsetCalibration.h.

Referenced by beginRun().

bool DTTTrigOffsetCalibration::doTTrigCorrection_
private

Definition at line 45 of file DTTTrigOffsetCalibration.h.

Referenced by beginRun(), and endJob().

TFile* DTTTrigOffsetCalibration::rootFile_
private
DTSegmentSelector DTTTrigOffsetCalibration::select_
private

Definition at line 42 of file DTTTrigOffsetCalibration.h.

Referenced by analyze().

std::string DTTTrigOffsetCalibration::theCalibChamber_
private

Definition at line 46 of file DTTTrigOffsetCalibration.h.

Referenced by analyze().

edm::InputTag DTTTrigOffsetCalibration::theRecHits4DLabel_
private

Definition at line 44 of file DTTTrigOffsetCalibration.h.

Referenced by analyze().

ChamberHistosMap DTTTrigOffsetCalibration::theT0SegHistoMap_
private

Definition at line 51 of file DTTTrigOffsetCalibration.h.

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

const DTTtrig* DTTTrigOffsetCalibration::tTrigMap_
private

Definition at line 50 of file DTTTrigOffsetCalibration.h.

Referenced by beginRun(), and endJob().