00001
00002
00003
00004
00005
00006
00007
00008 #include "CalibMuon/DTCalibration/plugins/DTTTrigCalibration.h"
00009 #include "CalibMuon/DTCalibration/interface/DTTimeBoxFitter.h"
00010 #include "CalibMuon/DTDigiSync/interface/DTTTrigSyncFactory.h"
00011 #include "CalibMuon/DTDigiSync/interface/DTTTrigBaseSync.h"
00012
00013 #include "FWCore/Framework/interface/EventSetup.h"
00014 #include "FWCore/Framework/interface/Event.h"
00015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00016 #include "FWCore/Framework/interface/ESHandle.h"
00017
00018 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00019 #include "DataFormats/MuonDetId/interface/DTWireId.h"
00020
00021 #include "CondFormats/DTObjects/interface/DTTtrig.h"
00022
00023 #include "DataFormats/GeometryVector/interface/GlobalPoint.h"
00024
00025 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00026
00027
00028 #include "CondFormats/DataRecord/interface/DTStatusFlagRcd.h"
00029 #include "CondFormats/DTObjects/interface/DTStatusFlag.h"
00030
00031
00032 #include "TFile.h"
00033 #include "TH1F.h"
00034 #include "TGraph.h"
00035
00036 class DTLayer;
00037
00038 using namespace std;
00039 using namespace edm;
00040
00041
00042
00043
00044
00045 DTTTrigCalibration::DTTTrigCalibration(const edm::ParameterSet& pset) {
00046
00047 debug = pset.getUntrackedParameter<bool>("debug");
00048
00049
00050 digiLabel = pset.getUntrackedParameter<string>("digiLabel");
00051
00052
00053 findTMeanAndSigma = pset.getUntrackedParameter<bool>("fitAndWrite", "false");
00054
00055
00056 maxTDCCounts = 5000 * pset.getUntrackedParameter<int>("tdcRescale", 1);
00057
00058 maxDigiPerLayer = pset.getUntrackedParameter<int>("maxDigiPerLayer", 10);
00059
00060
00061 string rootFileName = pset.getUntrackedParameter<string>("rootFileName");
00062 theFile = new TFile(rootFileName.c_str(), "RECREATE");
00063 theFile->cd();
00064 theFitter = new DTTimeBoxFitter();
00065 if(debug)
00066 theFitter->setVerbosity(1);
00067
00068 double sigmaFit = pset.getUntrackedParameter<double>("sigmaTTrigFit",10.);
00069 theFitter->setFitSigma(sigmaFit);
00070
00071 doSubtractT0 = pset.getUntrackedParameter<bool>("doSubtractT0","false");
00072
00073 if(doSubtractT0) {
00074 theSync = DTTTrigSyncFactory::get()->create(pset.getUntrackedParameter<string>("tTrigMode"),
00075 pset.getUntrackedParameter<ParameterSet>("tTrigModeConfig"));
00076 } else {
00077 theSync = 0;
00078 }
00079
00080 checkNoisyChannels = pset.getUntrackedParameter<bool>("checkNoisyChannels","false");
00081
00082 if(debug)
00083 cout << "[DTTTrigCalibration]Constructor called!" << endl;
00084 }
00085
00086
00087
00088
00089 DTTTrigCalibration::~DTTTrigCalibration(){
00090 if(debug)
00091 cout << "[DTTTrigCalibration]Destructor called!" << endl;
00092
00093
00094
00095
00096
00097
00098
00099
00100 theFile->Close();
00101 delete theFitter;
00102 }
00103
00104
00105
00107 void DTTTrigCalibration::analyze(const edm::Event & event, const edm::EventSetup& eventSetup) {
00108
00109 if(debug)
00110 cout << "[DTTTrigCalibration] #Event: " << event.id().event() << endl;
00111
00112
00113 Handle<DTDigiCollection> digis;
00114 event.getByLabel(digiLabel, digis);
00115
00116 ESHandle<DTStatusFlag> statusMap;
00117 if(checkNoisyChannels) {
00118
00119 eventSetup.get<DTStatusFlagRcd>().get(statusMap);
00120 }
00121
00122 if(doSubtractT0)
00123 theSync->setES(eventSetup);
00124
00125
00126 vector<DTChamberId> badChambers;
00127
00128
00129 DTDigiCollection::DigiRangeIterator dtLayerIt;
00130 for (dtLayerIt = digis->begin();
00131 dtLayerIt != digis->end();
00132 ++dtLayerIt){
00133
00134 const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;
00135
00136 const DTLayerId layerId = (*dtLayerIt).first;
00137 const DTSuperLayerId slId = layerId.superlayerId();
00138 const DTChamberId chId = slId.chamberId();
00139 bool badChamber=false;
00140
00141 if(debug)
00142 cout<<"----------- Layer "<<layerId<<" -------------"<<endl;
00143
00144
00145 for(vector<DTChamberId>::const_iterator chamber = badChambers.begin(); chamber != badChambers.end(); chamber++){
00146 if((*chamber) == chId){
00147 badChamber=true;
00148 break;
00149 }
00150 }
00151 if(badChamber) continue;
00152
00153
00154 if((digiRange.second - digiRange.first) > maxDigiPerLayer){
00155 if(debug)
00156 cout<<"Layer "<<layerId<<"has too many digis ("<<(digiRange.second - digiRange.first)<<")"<<endl;
00157 badChambers.push_back(chId);
00158 continue;
00159 }
00160
00161
00162 TH1F *hTBox = theHistoMap[slId];
00163 if(hTBox == 0) {
00164
00165 theFile->cd();
00166 hTBox = new TH1F(getTBoxName(slId).c_str(), "Time box (ns)", int(0.25*32.0*maxTDCCounts/25.0), 0, maxTDCCounts);
00167 if(debug)
00168 cout << " New Time Box: " << hTBox->GetName() << endl;
00169 theHistoMap[slId] = hTBox;
00170 }
00171 TH1F *hO = theOccupancyMap[layerId];
00172 if(hO == 0) {
00173
00174 theFile->cd();
00175 hO = new TH1F(getOccupancyName(layerId).c_str(), "Occupancy", 100, 0, 100);
00176 if(debug)
00177 cout << " New Time Box: " << hO->GetName() << endl;
00178 theOccupancyMap[layerId] = hO;
00179 }
00180
00181
00182 for (DTDigiCollection::const_iterator digi = digiRange.first;
00183 digi != digiRange.second;
00184 digi++) {
00185 const DTWireId wireId(layerId, (*digi).wire());
00186
00187
00188 if(checkNoisyChannels) {
00189 bool isNoisy = false;
00190 bool isFEMasked = false;
00191 bool isTDCMasked = false;
00192 bool isTrigMask = false;
00193 bool isDead = false;
00194 bool isNohv = false;
00195 statusMap->cellStatus(wireId, isNoisy, isFEMasked, isTDCMasked, isTrigMask, isDead, isNohv);
00196 if(isNoisy) {
00197 if(debug)
00198 cout << "Wire: " << wireId << " is noisy, skipping!" << endl;
00199 continue;
00200 }
00201 }
00202 theFile->cd();
00203 double offset = 0;
00204 if(doSubtractT0) {
00205 const DTLayer* layer = 0;
00206 const GlobalPoint glPt;
00207 offset = theSync->offset(layer, wireId, glPt);
00208 }
00209 hTBox->Fill((*digi).time()-offset);
00210 if(debug) {
00211 cout << " Filling Time Box: " << hTBox->GetName() << endl;
00212 cout << " offset (ns): " << offset << endl;
00213 cout << " time(ns): " << (*digi).time()-offset<< endl;
00214 }
00215 hO->Fill((*digi).wire());
00216 }
00217 }
00218 }
00219
00220
00221 void DTTTrigCalibration::endJob() {
00222 if(debug)
00223 cout << "[DTTTrigCalibration]Writing histos to file!" << endl;
00224
00225
00226 theFile->cd();
00227 for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
00228 slHisto != theHistoMap.end();
00229 slHisto++) {
00230 (*slHisto).second->Write();
00231 }
00232 for(map<DTLayerId, TH1F*>::const_iterator slHisto = theOccupancyMap.begin();
00233 slHisto != theOccupancyMap.end();
00234 slHisto++) {
00235 (*slHisto).second->Write();
00236 }
00237
00238 if(findTMeanAndSigma) {
00239
00240 DTTtrig* tTrig = new DTTtrig();
00241
00242
00243 for(map<DTSuperLayerId, TH1F*>::const_iterator slHisto = theHistoMap.begin();
00244 slHisto != theHistoMap.end();
00245 slHisto++) {
00246 pair<double, double> meanAndSigma = theFitter->fitTimeBox((*slHisto).second);
00247 tTrig->set((*slHisto).first,
00248 meanAndSigma.first,
00249 meanAndSigma.second,
00250 DTTimeUnits::ns);
00251
00252 if(debug) {
00253 cout << " SL: " << (*slHisto).first
00254 << " mean = " << meanAndSigma.first
00255 << " sigma = " << meanAndSigma.second << endl;
00256 }
00257 }
00258
00259
00260 dumpTTrigMap(tTrig);
00261
00262
00263 plotTTrig(tTrig);
00264
00265 if(debug)
00266 cout << "[DTTTrigCalibration]Writing ttrig object to DB!" << endl;
00267
00268
00269
00270 string tTrigRecord = "DTTtrigRcd";
00271
00272
00273 DTCalibDBUtils::writeToDB(tTrigRecord, tTrig);
00274 }
00275
00276 }
00277
00278
00279
00280
00281 string DTTTrigCalibration::getTBoxName(const DTSuperLayerId& slId) const {
00282 string histoName;
00283 stringstream theStream;
00284 theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
00285 << "_SL" << slId.superlayer() << "_hTimeBox";
00286 theStream >> histoName;
00287 return histoName;
00288 }
00289
00290 string DTTTrigCalibration::getOccupancyName(const DTLayerId& slId) const {
00291 string histoName;
00292 stringstream theStream;
00293 theStream << "Ch_" << slId.wheel() << "_" << slId.station() << "_" << slId.sector()
00294 << "_SL" << slId.superlayer() << "_L"<< slId.layer() <<"_Occupancy";
00295 theStream >> histoName;
00296 return histoName;
00297 }
00298
00299
00300 void DTTTrigCalibration::dumpTTrigMap(const DTTtrig* tTrig) const {
00301 static const double convToNs = 25./32.;
00302 for(DTTtrig::const_iterator ttrig = tTrig->begin();
00303 ttrig != tTrig->end(); ttrig++) {
00304 cout << "Wh: " << (*ttrig).first.wheelId
00305 << " St: " << (*ttrig).first.stationId
00306 << " Sc: " << (*ttrig).first.sectorId
00307 << " Sl: " << (*ttrig).first.slId
00308 << " TTrig mean (ns): " << (*ttrig).second.tTrig * convToNs
00309 << " TTrig sigma (ns): " << (*ttrig).second.tTrms * convToNs<< endl;
00310 }
00311 }
00312
00313
00314 void DTTTrigCalibration::plotTTrig(const DTTtrig* tTrig) const {
00315
00316 TH1F* tTrig_YB1_Se10 = new TH1F("tTrig_YB1_Se10","tTrig YB1_Se10",15,1,16);
00317 TH1F* tTrig_YB2_Se10 = new TH1F("tTrig_YB2_Se10","tTrig YB2_Se10",15,1,16);
00318 TH1F* tTrig_YB2_Se11 = new TH1F("tTrig_YB2_Se11","tTrig YB2_Se11",12,1,13);
00319
00320 static const double convToNs = 25./32.;
00321 for(DTTtrig::const_iterator ttrig = tTrig->begin();
00322 ttrig != tTrig->end(); ttrig++) {
00323
00324
00325 float tTrigValue=0;
00326 float tTrmsValue=0;
00327 if ((*ttrig).second.tTrig * convToNs > 0 &&
00328 (*ttrig).second.tTrig * convToNs < 32000 ) {
00329 tTrigValue = (*ttrig).second.tTrig * convToNs;
00330 tTrmsValue = (*ttrig).second.tTrms * convToNs;
00331 }
00332
00333 int binx;
00334 string binLabel;
00335 stringstream binLabelStream;
00336 if ((*ttrig).first.sectorId != 14) {
00337 binx = ((*ttrig).first.stationId-1)*3 + (*ttrig).first.slId;
00338 binLabelStream << "MB"<<(*ttrig).first.stationId<<"_SL"<<(*ttrig).first.slId;
00339 }
00340 else {
00341 binx = 12 + (*ttrig).first.slId;
00342 binLabelStream << "MB14_SL"<<(*ttrig).first.slId;
00343 }
00344 binLabelStream >> binLabel;
00345
00346 if ((*ttrig).first.wheelId == 2) {
00347 if ((*ttrig).first.sectorId == 10 || (*ttrig).first.sectorId == 14) {
00348 tTrig_YB2_Se10->Fill( binx,tTrigValue);
00349 tTrig_YB2_Se10->SetBinError( binx, tTrmsValue);
00350 tTrig_YB2_Se10->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
00351 tTrig_YB2_Se10->GetYaxis()->SetTitle("ns");
00352 }
00353 else {
00354 tTrig_YB2_Se11->Fill( binx,tTrigValue);
00355 tTrig_YB2_Se11->SetBinError( binx,tTrmsValue);
00356 tTrig_YB2_Se11->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
00357 tTrig_YB2_Se11->GetYaxis()->SetTitle("ns");
00358 }
00359 }
00360 else {
00361 tTrig_YB1_Se10->Fill( binx,tTrigValue);
00362 tTrig_YB1_Se10->SetBinError( binx,tTrmsValue);
00363 tTrig_YB1_Se10->GetXaxis()->SetBinLabel(binx,binLabel.c_str());
00364 tTrig_YB1_Se10->GetYaxis()->SetTitle("ns");
00365 }
00366 }
00367
00368 tTrig_YB1_Se10->Write();
00369 tTrig_YB2_Se10->Write();
00370 tTrig_YB2_Se11->Write();
00371
00372 }