Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:34:57

0001 
0002 /*
0003  *  See header file for a description of this class.
0004  *
0005  *  \author A. Vilela Pereira
0006  */
0007 
0008 #include "DTTTrigOffsetCalibration.h"
0009 
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Framework/interface/EventSetup.h"
0012 #include "FWCore/Framework/interface/Event.h"
0013 #include "FWCore/Framework/interface/ESHandle.h"
0014 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0015 #include "FWCore/Utilities/interface/InputTag.h"
0016 
0017 #include "Geometry/DTGeometry/interface/DTGeometry.h"
0018 #include "Geometry/Records/interface/MuonGeometryRecord.h"
0019 
0020 #include "CondFormats/DTObjects/interface/DTTtrig.h"
0021 #include "CondFormats/DataRecord/interface/DTTtrigRcd.h"
0022 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
0023 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
0024 
0025 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
0026 #include "CalibMuon/DTCalibration/interface/DTSegmentSelector.h"
0027 
0028 #include <string>
0029 #include <sstream>
0030 #include "TFile.h"
0031 #include "TH1F.h"
0032 
0033 using namespace std;
0034 using namespace edm;
0035 
0036 DTTTrigOffsetCalibration::DTTTrigOffsetCalibration(const ParameterSet& pset)
0037     : theRecHits4DToken_(consumes<DTRecSegment4DCollection>(pset.getParameter<InputTag>("recHits4DLabel"))),
0038       doTTrigCorrection_(pset.getUntrackedParameter<bool>("doT0SegCorrection", false)),
0039       theCalibChamber_(pset.getUntrackedParameter<string>("calibChamber", "All")),
0040       ttrigToken_(
0041           esConsumes<edm::Transition::BeginRun>(edm::ESInputTag("", pset.getUntrackedParameter<string>("dbLabel")))),
0042       dtGeomToken_(esConsumes()) {
0043   LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Constructor called!";
0044 
0045   edm::ConsumesCollector collector(consumesCollector());
0046   select_ = new DTSegmentSelector(pset, collector);
0047 
0048   // the root file which will contain the histos
0049   string rootFileName = pset.getUntrackedParameter<string>("rootFileName", "DTT0SegHistos.root");
0050   rootFile_ = new TFile(rootFileName.c_str(), "RECREATE");
0051   rootFile_->cd();
0052 }
0053 
0054 void DTTTrigOffsetCalibration::beginRun(const edm::Run& run, const edm::EventSetup& setup) {
0055   if (doTTrigCorrection_) {
0056     ESHandle<DTTtrig> tTrig;
0057     tTrig = setup.getHandle(ttrigToken_);
0058     tTrigMap_ = &setup.getData(ttrigToken_);
0059     LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration]: TTrig version: " << tTrig->version() << endl;
0060   }
0061 }
0062 
0063 DTTTrigOffsetCalibration::~DTTTrigOffsetCalibration() {
0064   rootFile_->Close();
0065   LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Destructor called!";
0066 }
0067 
0068 void DTTTrigOffsetCalibration::analyze(const Event& event, const EventSetup& eventSetup) {
0069   rootFile_->cd();
0070   DTChamberId chosenChamberId;
0071 
0072   if (theCalibChamber_ != "All") {
0073     stringstream linestr;
0074     int selWheel, selStation, selSector;
0075     linestr << theCalibChamber_;
0076     linestr >> selWheel >> selStation >> selSector;
0077     chosenChamberId = DTChamberId(selWheel, selStation, selSector);
0078     LogVerbatim("Calibration") << " chosen chamber " << chosenChamberId << endl;
0079   }
0080 
0081   // Get the DT Geometry
0082   ESHandle<DTGeometry> dtGeom;
0083   dtGeom = eventSetup.getHandle(dtGeomToken_);
0084 
0085   // Get the rechit collection from the event
0086   const Handle<DTRecSegment4DCollection>& all4DSegments = event.getHandle(theRecHits4DToken_);
0087 
0088   // Loop over segments by chamber
0089   DTRecSegment4DCollection::id_iterator chamberIdIt;
0090   for (chamberIdIt = all4DSegments->id_begin(); chamberIdIt != all4DSegments->id_end(); ++chamberIdIt) {
0091     // Get the chamber from the setup
0092     const DTChamber* chamber = dtGeom->chamber(*chamberIdIt);
0093     LogTrace("Calibration") << "Chamber Id: " << *chamberIdIt;
0094 
0095     // Book histos
0096     if (theT0SegHistoMap_.find(*chamberIdIt) == theT0SegHistoMap_.end()) {
0097       bookHistos(*chamberIdIt);
0098     }
0099 
0100     // Calibrate just the chosen chamber/s
0101     if ((theCalibChamber_ != "All") && ((*chamberIdIt) != chosenChamberId))
0102       continue;
0103 
0104     // Get the range for the corresponding ChamberId
0105     DTRecSegment4DCollection::range range = all4DSegments->get((*chamberIdIt));
0106 
0107     // Loop over the rechits of this DetUnit
0108     for (DTRecSegment4DCollection::const_iterator segment = range.first; segment != range.second; ++segment) {
0109       LogTrace("Calibration") << "Segment local pos (in chamber RF): " << (*segment).localPosition()
0110                               << "\nSegment global pos: " << chamber->toGlobal((*segment).localPosition());
0111 
0112       if (!(*select_)(*segment, event, eventSetup))
0113         continue;
0114 
0115       // Fill t0-seg values
0116       if ((*segment).hasPhi()) {
0117         //if( segment->phiSegment()->ist0Valid() ){
0118         if ((segment->phiSegment()->t0()) != 0.00) {
0119           (theT0SegHistoMap_[*chamberIdIt])[0]->Fill(segment->phiSegment()->t0());
0120         }
0121       }
0122       if ((*segment).hasZed()) {
0123         //if( segment->zSegment()->ist0Valid() ){
0124         if ((segment->zSegment()->t0()) != 0.00) {
0125           (theT0SegHistoMap_[*chamberIdIt])[1]->Fill(segment->zSegment()->t0());
0126         }
0127       }
0128     }  // DTRecSegment4DCollection::const_iterator segment
0129   }  // DTRecSegment4DCollection::id_iterator chamberIdIt
0130 }  // DTTTrigOffsetCalibration::analyze
0131 
0132 void DTTTrigOffsetCalibration::endJob() {
0133   rootFile_->cd();
0134 
0135   LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Writing histos to file!" << endl;
0136 
0137   for (ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end();
0138        ++itChHistos) {
0139     for (vector<TH1F*>::const_iterator itHist = (*itChHistos).second.begin(); itHist != (*itChHistos).second.end();
0140          ++itHist)
0141       (*itHist)->Write();
0142   }
0143 
0144   if (doTTrigCorrection_) {
0145     // Create the object to be written to DB
0146     DTTtrig* tTrig = new DTTtrig();
0147 
0148     for (ChamberHistosMap::const_iterator itChHistos = theT0SegHistoMap_.begin(); itChHistos != theT0SegHistoMap_.end();
0149          ++itChHistos) {
0150       DTChamberId chId = itChHistos->first;
0151       // Get SuperLayerId's for each ChamberId
0152       vector<DTSuperLayerId> slIds;
0153       slIds.push_back(DTSuperLayerId(chId, 1));
0154       slIds.push_back(DTSuperLayerId(chId, 3));
0155       if (chId.station() != 4)
0156         slIds.push_back(DTSuperLayerId(chId, 2));
0157 
0158       for (vector<DTSuperLayerId>::const_iterator itSl = slIds.begin(); itSl != slIds.end(); ++itSl) {
0159         // Get old values from DB
0160         float ttrigMean = 0;
0161         float ttrigSigma = 0;
0162         float kFactor = 0;
0163         tTrigMap_->get(*itSl, ttrigMean, ttrigSigma, kFactor, DTTimeUnits::ns);
0164         //FIXME: verify if values make sense
0165         // Set new values
0166         float ttrigMeanNew = ttrigMean;
0167         float ttrigSigmaNew = ttrigSigma;
0168         float t0SegMean =
0169             (itSl->superLayer() != 2) ? itChHistos->second[0]->GetMean() : itChHistos->second[1]->GetMean();
0170 
0171         float kFactorNew = (kFactor * ttrigSigma + t0SegMean) / ttrigSigma;
0172 
0173         tTrig->set(*itSl, ttrigMeanNew, ttrigSigmaNew, kFactorNew, DTTimeUnits::ns);
0174       }
0175     }
0176     LogVerbatim("Calibration") << "[DTTTrigOffsetCalibration] Writing ttrig object to DB!" << endl;
0177     // Write the object to DB
0178     string tTrigRecord = "DTTtrigRcd";
0179     DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
0180   }
0181 }
0182 
0183 // Book a set of histograms for a given Chamber
0184 void DTTTrigOffsetCalibration::bookHistos(DTChamberId chId) {
0185   LogTrace("Calibration") << "   Booking histos for Chamber: " << chId;
0186 
0187   // Compose the chamber name
0188   std::string wheel = std::to_string(chId.wheel());
0189   std::string station = std::to_string(chId.station());
0190   std::string sector = std::to_string(chId.sector());
0191 
0192   string chHistoName = "_W" + wheel + "_St" + station + "_Sec" + sector;
0193 
0194   vector<TH1F*> histos;
0195   // Note the order matters
0196   histos.push_back(new TH1F(("hRPhiSegT0" + chHistoName).c_str(), "t0 from Phi segments", 500, -60., 60.));
0197   if (chId.station() != 4)
0198     histos.push_back(new TH1F(("hRZSegT0" + chHistoName).c_str(), "t0 from Z segments", 500, -60., 60.));
0199 
0200   theT0SegHistoMap_[chId] = histos;
0201 }