CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/CalibMuon/DTCalibration/plugins/DTT0Calibration.cc

Go to the documentation of this file.
00001 /*
00002  *  See header file for a description of this class.
00003  *
00004  *  $Date: 2012/02/17 16:39:06 $
00005  *  $Revision: 1.11.2.1 $
00006  *  \author S. Bolognesi - INFN Torino
00007  */
00008 #include "CalibMuon/DTCalibration/plugins/DTT0Calibration.h"
00009 #include "CalibMuon/DTCalibration/interface/DTCalibDBUtils.h"
00010 
00011 #include "FWCore/Framework/interface/Event.h"
00012 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00013 
00014 #include "Geometry/Records/interface/MuonGeometryRecord.h"
00015 #include "Geometry/Records/interface/MuonNumberingRecord.h"
00016 
00017 #include "DataFormats/DTDigi/interface/DTDigiCollection.h"
00018 #include "CondFormats/DTObjects/interface/DTT0.h"
00019 
00020 #include "TH1I.h"
00021 #include "TFile.h"
00022 #include "TKey.h"
00023 
00024 using namespace std;
00025 using namespace edm;
00026 // using namespace cond;
00027 
00028 // Constructor
00029 DTT0Calibration::DTT0Calibration(const edm::ParameterSet& pset) {
00030   // Get the debug parameter for verbose output
00031   debug = pset.getUntrackedParameter<bool>("debug");
00032   if(debug) 
00033     cout << "[DTT0Calibration]Constructor called!" << endl;
00034 
00035   // Get the label to retrieve digis from the event
00036   digiLabel = pset.getUntrackedParameter<string>("digiLabel");
00037 
00038   // The root file which contain the histos per layer
00039   string rootFileName = pset.getUntrackedParameter<string>("rootFileName","DTT0PerLayer.root");
00040   theFile = new TFile(rootFileName.c_str(), "RECREATE");
00041  
00042   theCalibWheel =  pset.getUntrackedParameter<string>("calibWheel", "All"); //FIXME amke a vector of integer instead of a string
00043   if(theCalibWheel != "All") {
00044     stringstream linestr;
00045     int selWheel;
00046     linestr << theCalibWheel;
00047     linestr >> selWheel;
00048     cout << "[DTT0CalibrationPerLayer] chosen wheel " << selWheel << endl;
00049   }
00050 
00051   // Sector/s to calibrate
00052   theCalibSector =  pset.getUntrackedParameter<string>("calibSector", "All"); //FIXME amke a vector of integer instead of a string
00053   if(theCalibSector != "All") {
00054     stringstream linestr;
00055     int selSector;
00056     linestr << theCalibSector;
00057     linestr >> selSector;
00058     cout << "[DTT0CalibrationPerLayer] chosen sector " << selSector << endl;
00059   }
00060 
00061   vector<string> defaultCell;
00062   defaultCell.push_back("None");
00063   cellsWithHistos = pset.getUntrackedParameter<vector<string> >("cellsWithHisto", defaultCell);
00064   for(vector<string>::const_iterator cell = cellsWithHistos.begin(); cell != cellsWithHistos.end(); cell++){
00065     if((*cell) != "None"){
00066       stringstream linestr;
00067       int wheel,sector,station,sl,layer,wire;
00068       linestr << (*cell);
00069       linestr >> wheel >> sector >> station >> sl >> layer >> wire;
00070       wireIdWithHistos.push_back(DTWireId(wheel,station,sector,sl,layer,wire));
00071     }
00072   }
00073 
00074   hT0SectorHisto=0;
00075 
00076   nevents=0;
00077   eventsForLayerT0 = pset.getParameter<unsigned int>("eventsForLayerT0");
00078   eventsForWireT0 = pset.getParameter<unsigned int>("eventsForWireT0");
00079   rejectDigiFromPeak = pset.getParameter<unsigned int>("rejectDigiFromPeak");
00080   tpPeakWidth = pset.getParameter<double>("tpPeakWidth");
00081 }
00082 
00083 // Destructor
00084 DTT0Calibration::~DTT0Calibration(){
00085   if(debug) 
00086     cout << "[DTT0Calibration]Destructor called!" << endl;
00087 
00088   theFile->Close();
00089 }
00090 
00092 void DTT0Calibration::analyze(const edm::Event & event, const edm::EventSetup& eventSetup) {
00093   if(debug || event.id().event() % 500==0)
00094     cout << "--- [DTT0Calibration] Analysing Event: #Run: " << event.id().run()
00095          << " #Event: " << event.id().event() << endl;
00096   nevents++;
00097 
00098   // Get the digis from the event
00099   Handle<DTDigiCollection> digis; 
00100   event.getByLabel(digiLabel, digis);
00101 
00102   // Get the DT Geometry
00103   eventSetup.get<MuonGeometryRecord>().get(dtGeom);
00104 
00105   // Iterate through all digi collections ordered by LayerId   
00106   DTDigiCollection::DigiRangeIterator dtLayerIt;
00107   for (dtLayerIt = digis->begin();
00108        dtLayerIt != digis->end();
00109        ++dtLayerIt){
00110     // Get the iterators over the digis associated with this LayerId
00111     const DTDigiCollection::Range& digiRange = (*dtLayerIt).second;
00112   
00113     // Get the layerId
00114     const DTLayerId layerId = (*dtLayerIt).first; //FIXME: check to be in the right sector
00115 
00116     if((theCalibWheel != "All") && (layerId.superlayerId().chamberId().wheel() != selWheel))
00117       continue;
00118     if((theCalibSector != "All") && (layerId.superlayerId().chamberId().sector() != selSector))
00119       continue;
00120  
00121     //if(debug) {
00122     //  cout << "Layer " << layerId<<" with "<<distance(digiRange.first, digiRange.second)<<" digi"<<endl;
00123     //}
00124 
00125     // Loop over all digis in the given layer
00126     for (DTDigiCollection::const_iterator digi = digiRange.first;
00127          digi != digiRange.second;
00128          digi++) {
00129       double t0 = (*digi).countsTDC();
00130 
00131       //Use first bunch of events to fill t0 per layer
00132       if(nevents < eventsForLayerT0){
00133         //Get the per-layer histo from the map
00134         TH1I *hT0LayerHisto = theHistoLayerMap[layerId];
00135         //If it doesn't exist, book it
00136         if(hT0LayerHisto == 0){
00137           theFile->cd();
00138           hT0LayerHisto = new TH1I(getHistoName(layerId).c_str(),
00139                                    "T0 from pulses by layer (TDC counts, 1 TDC count = 0.781 ns)",
00140                                    200, t0-100, t0+100);
00141           if(debug)
00142             cout << "  New T0 per Layer Histo: " << hT0LayerHisto->GetName() << endl;
00143           theHistoLayerMap[layerId] = hT0LayerHisto;
00144         }
00145     
00146         //Fill the histos
00147         theFile->cd();
00148         if(hT0LayerHisto != 0) {
00149           //  if(debug)
00150           // cout<<"Filling histo "<<hT0LayerHisto->GetName()<<" with digi "<<t0<<" TDC counts"<<endl;
00151           hT0LayerHisto->Fill(t0);
00152         }
00153       }
00154 
00155       //Use all the remaining events to compute t0 per wire
00156       if(nevents>eventsForLayerT0){
00157         // Get the wireId
00158         const DTWireId wireId(layerId, (*digi).wire());
00159         if(debug) {
00160           cout << "   Wire: " << wireId << endl
00161                << "       time (TDC counts): " << (*digi).countsTDC()<< endl;
00162         }   
00163 
00164         //Fill the histos per wire for the chosen cells
00165         vector<DTWireId>::iterator it_wire = find(wireIdWithHistos.begin(),wireIdWithHistos.end(),wireId);
00166         if(it_wire != wireIdWithHistos.end()){
00167           if(theHistoWireMap.find(wireId) == theHistoWireMap.end()){
00168             theHistoWireMap[wireId] = new TH1I(getHistoName(wireId).c_str(),"T0 from pulses by wire (TDC counts, 1 TDC count = 0.781 ns)",7000,0,7000);
00169             if(debug) cout << "  New T0 per wire Histo: " << (theHistoWireMap[wireId])->GetName() << endl;
00170           }
00171           if(theHistoWireMap_ref.find(wireId) == theHistoWireMap_ref.end()){
00172             theHistoWireMap_ref[wireId] = new TH1I((getHistoName(wireId) + "_ref").c_str(),"T0 from pulses by wire (TDC counts, 1 TDC count = 0.781 ns)",7000,0,7000);
00173             if(debug) cout << "  New T0 per wire Histo: " << (theHistoWireMap_ref[wireId])->GetName() << endl;
00174           }
00175 
00176           TH1I* hT0WireHisto = theHistoWireMap[wireId];
00177           //Fill the histos
00178           theFile->cd();
00179           if(hT0WireHisto) hT0WireHisto->Fill(t0);
00180         }
00181 
00182         //Check the tzero has reasonable value
00183         if(abs(hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin()) - t0) > rejectDigiFromPeak){
00184           if(debug)
00185             cout<<"digi skipped because t0 per sector "<<hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin())<<endl;
00186           continue;
00187         }
00188 
00189         //Use second bunch of events to compute a t0 reference per wire
00190         if(nevents< (eventsForLayerT0 + eventsForWireT0)){
00191           //Fill reference wire histos
00192           if(it_wire != wireIdWithHistos.end()){
00193             TH1I* hT0WireHisto_ref = theHistoWireMap_ref[wireId];
00194             theFile->cd();
00195             if(hT0WireHisto_ref) hT0WireHisto_ref->Fill(t0); 
00196           } 
00197           if(!nDigiPerWire_ref[wireId]){
00198             mK_ref[wireId] = 0;
00199           }
00200           nDigiPerWire_ref[wireId] = nDigiPerWire_ref[wireId] + 1;
00201           mK_ref[wireId] = mK_ref[wireId] + (t0-mK_ref[wireId])/nDigiPerWire_ref[wireId];
00202         }
00203         //Use last all the remaining events to compute the mean and sigma t0 per wire
00204         else if(nevents>(eventsForLayerT0 + eventsForWireT0)){
00205           if(abs(t0-mK_ref[wireId]) > tpPeakWidth)
00206             continue;
00207           if(!nDigiPerWire[wireId]){
00208             theAbsoluteT0PerWire[wireId] = 0;
00209             qK[wireId] = 0;
00210             mK[wireId] = 0;
00211           }
00212           nDigiPerWire[wireId] = nDigiPerWire[wireId] + 1;
00213           theAbsoluteT0PerWire[wireId] = theAbsoluteT0PerWire[wireId] + t0;
00214           //theSigmaT0PerWire[wireId] = theSigmaT0PerWire[wireId] + (t0*t0);
00215           qK[wireId] = qK[wireId] + ((nDigiPerWire[wireId]-1)*(t0-mK[wireId])*(t0-mK[wireId])/nDigiPerWire[wireId]);
00216           mK[wireId] = mK[wireId] + (t0-mK[wireId])/nDigiPerWire[wireId];
00217         }
00218       }//end if(nevents>1000)
00219     }//end loop on digi
00220   }//end loop on layer
00221 
00222   //Use the t0 per layer histos to have an indication about the t0 position 
00223   if(nevents == eventsForLayerT0){
00224     for(map<DTLayerId, TH1I*>::const_iterator lHisto = theHistoLayerMap.begin();
00225         lHisto != theHistoLayerMap.end();
00226         lHisto++){
00227       if(debug)
00228         cout<<"Reading histogram "<<(*lHisto).second->GetName()<<" with mean "<<(*lHisto).second->GetMean()<<" and RMS "<<(*lHisto).second->GetRMS();
00229 
00230       //Take the mean as a first t0 estimation
00231       if((*lHisto).second->GetRMS()<5.0){
00232         if(hT0SectorHisto == 0){
00233           hT0SectorHisto = new TH1D("hT0AllLayerOfSector","T0 from pulses per layer in sector", 
00234                                     //20, (*lHisto).second->GetMean()-100, (*lHisto).second->GetMean()+100);
00235                                     700, 0, 7000);
00236         }
00237         if(debug)
00238           cout<<" accepted"<<endl;
00239         hT0SectorHisto->Fill((*lHisto).second->GetMean());
00240       }
00241       //Take the mean of noise + 400ns as a first t0 estimation
00242       // if((*lHisto).second->GetRMS()>10.0 && ((*lHisto).second->GetRMS()<15.0)){
00243 //      double t0_estim = (*lHisto).second->GetMean() + 400;
00244 //      if(hT0SectorHisto == 0){
00245 //        hT0SectorHisto = new TH1D("hT0AllLayerOfSector","T0 from pulses per layer in sector", 
00246 //                                  //20, t0_estim-100, t0_estim+100);
00247 //                                  700, 0, 7000);
00248 //      }
00249 //      if(debug)
00250 //        cout<<" accepted + 400ns"<<endl;
00251 //      hT0SectorHisto->Fill((*lHisto).second->GetMean() + 400);
00252 //       }
00253       if(debug)
00254         cout<<endl;
00255 
00256       theT0LayerMap[(*lHisto).second->GetName()] = (*lHisto).second->GetMean();
00257       theSigmaT0LayerMap[(*lHisto).second->GetName()] = (*lHisto).second->GetRMS();
00258     }
00259     if(!hT0SectorHisto){
00260       cout<<"[DTT0Calibration]: All the t0 per layer are still uncorrect: trying with greater number of events"<<endl;
00261       eventsForLayerT0 = eventsForLayerT0*2;
00262       return;
00263     }
00264     if(debug)
00265       cout<<"[DTT0Calibration] t0 reference for this sector "<<
00266         hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin())<<endl;
00267   } 
00268 }
00269 
00270 
00271 void DTT0Calibration::endJob() {
00272 
00273   DTT0* t0s = new DTT0();
00274   DTT0* t0sWRTChamber = new DTT0();
00275 
00276   if(debug) 
00277     cout << "[DTT0CalibrationPerLayer]Writing histos to file!" << endl;
00278 
00279   theFile->cd();
00280   hT0SectorHisto->Write();
00281   for(map<DTWireId, TH1I*>::const_iterator wHisto = theHistoWireMap.begin();
00282       wHisto != theHistoWireMap.end();
00283       wHisto++) {
00284     (*wHisto).second->Write(); 
00285   }
00286   for(map<DTWireId, TH1I*>::const_iterator wHisto = theHistoWireMap_ref.begin();
00287       wHisto != theHistoWireMap_ref.end();
00288       wHisto++) {
00289     (*wHisto).second->Write();
00290   }
00291   for(map<DTLayerId, TH1I*>::const_iterator lHisto = theHistoLayerMap.begin();
00292       lHisto != theHistoLayerMap.end();
00293       lHisto++) {
00294     (*lHisto).second->Write(); 
00295   }  
00296 
00297   if(debug) 
00298     cout << "[DTT0Calibration] Compute and store t0 and sigma per wire" << endl;
00299 
00300   for(map<DTWireId, double>::const_iterator wiret0 = theAbsoluteT0PerWire.begin();
00301       wiret0 != theAbsoluteT0PerWire.end();
00302       wiret0++){
00303     if(nDigiPerWire[(*wiret0).first]){
00304       double t0 = (*wiret0).second/nDigiPerWire[(*wiret0).first];
00305       theRelativeT0PerWire[(*wiret0).first] = t0 - hT0SectorHisto->GetBinCenter(hT0SectorHisto->GetMaximumBin());
00306       cout<<"Wire "<<(*wiret0).first<<" has    t0 "<<t0<<"(absolute) "<<theRelativeT0PerWire[(*wiret0).first]<<"(relative)";
00307 
00308       //theSigmaT0PerWire[(*wiret0).first] = sqrt((theSigmaT0PerWire[(*wiret0).first] / nDigiPerWire[(*wiret0).first]) - t0*t0);
00309       theSigmaT0PerWire[(*wiret0).first] = sqrt(qK[(*wiret0).first]/nDigiPerWire[(*wiret0).first]);
00310       cout<<"    sigma "<<theSigmaT0PerWire[(*wiret0).first]<<endl;
00311     }
00312     else{
00313       cout<<"[DTT0Calibration] ERROR: no digis in wire "<<(*wiret0).first<<endl;
00314       abort();
00315     }
00316   }
00317 
00319   // Get all the sls from the setup
00320   const vector<DTSuperLayer*> superLayers = dtGeom->superLayers();     
00321   // Loop over all SLs
00322   for(vector<DTSuperLayer*>::const_iterator  sl = superLayers.begin();
00323       sl != superLayers.end(); sl++) {
00324 
00325   
00326     //Compute mean for odd and even superlayers
00327     double oddLayersMean=0;
00328     double evenLayersMean=0; 
00329     double oddLayersDen=0;
00330     double evenLayersDen=0;
00331     for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00332         wiret0 != theRelativeT0PerWire.end();
00333         wiret0++){
00334       if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00335         if(debug)
00336           cout<<"[DTT0Calibration] Superlayer "<<(*sl)->id()
00337               <<"layer " <<(*wiret0).first.layerId().layer()<<" with "<<(*wiret0).second<<endl;
00338         if(((*wiret0).first.layerId().layer()) % 2){
00339           oddLayersMean = oddLayersMean + (*wiret0).second;
00340           oddLayersDen++;
00341         }
00342         else{
00343           evenLayersMean = evenLayersMean + (*wiret0).second;
00344           evenLayersDen++;
00345         }
00346       }
00347     }
00348     oddLayersMean = oddLayersMean/oddLayersDen;
00349     evenLayersMean = evenLayersMean/evenLayersDen;
00350     if(debug && oddLayersMean)
00351       cout<<"[DTT0Calibration] Relative T0 mean for  odd layers "<<oddLayersMean<<"  even layers"<<evenLayersMean<<endl;
00352 
00353     //Compute sigma for odd and even superlayers
00354     double oddLayersSigma=0;
00355     double evenLayersSigma=0;
00356     for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00357         wiret0 != theRelativeT0PerWire.end();
00358         wiret0++){
00359       if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00360         if(((*wiret0).first.layerId().layer()) % 2){
00361           oddLayersSigma = oddLayersSigma + ((*wiret0).second - oddLayersMean) * ((*wiret0).second - oddLayersMean);
00362         }
00363         else{
00364           evenLayersSigma = evenLayersSigma + ((*wiret0).second - evenLayersMean) * ((*wiret0).second - evenLayersMean);
00365         }
00366       }
00367     }
00368     oddLayersSigma = oddLayersSigma/oddLayersDen;
00369     evenLayersSigma = evenLayersSigma/evenLayersDen;
00370     oddLayersSigma = sqrt(oddLayersSigma);
00371     evenLayersSigma = sqrt(evenLayersSigma);
00372 
00373     if(debug && oddLayersMean)
00374       cout<<"[DTT0Calibration] Relative T0 sigma for  odd layers "<<oddLayersSigma<<"  even layers"<<evenLayersSigma<<endl;
00375 
00376     //Recompute the mean for odd and even superlayers discarding fluctations
00377     double oddLayersFinalMean=0; 
00378     double evenLayersFinalMean=0;
00379     for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00380         wiret0 != theRelativeT0PerWire.end();
00381         wiret0++){
00382       if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00383         if(((*wiret0).first.layerId().layer()) % 2){
00384           if(abs((*wiret0).second - oddLayersMean) < (2*oddLayersSigma))
00385             oddLayersFinalMean = oddLayersFinalMean + (*wiret0).second;
00386         }
00387         else{
00388           if(abs((*wiret0).second - evenLayersMean) < (2*evenLayersSigma))
00389             evenLayersFinalMean = evenLayersFinalMean + (*wiret0).second;
00390         }
00391       }
00392     }
00393     oddLayersFinalMean = oddLayersFinalMean/oddLayersDen;
00394     evenLayersFinalMean = evenLayersFinalMean/evenLayersDen;
00395     if(debug && oddLayersMean)
00396       cout<<"[DTT0Calibration] Final relative T0 mean for  odd layers "<<oddLayersFinalMean<<"  even layers"<<evenLayersFinalMean<<endl;
00397 
00398     for(map<DTWireId, double>::const_iterator wiret0 = theRelativeT0PerWire.begin();
00399         wiret0 != theRelativeT0PerWire.end();
00400         wiret0++){
00401       if((*wiret0).first.layerId().superlayerId() == (*sl)->id()){
00402         double t0=-999;
00403         if(((*wiret0).first.layerId().layer()) % 2)
00404           t0 = (*wiret0).second + (evenLayersFinalMean - oddLayersFinalMean);
00405         else
00406           t0 = (*wiret0).second;
00407 
00408         cout<<"[DTT0Calibration] Wire "<<(*wiret0).first<<" has    t0 "<<(*wiret0).second<<" (relative, after even-odd layer corrections)  "
00409             <<"    sigma "<<theSigmaT0PerWire[(*wiret0).first]<<endl;
00410         //Store the results into DB
00411         t0s->set((*wiret0).first, t0, theSigmaT0PerWire[(*wiret0).first],DTTimeUnits::counts); 
00412       }
00413     }
00414   }
00415   
00417   if(debug) 
00418     cout << "[DTT0Calibration]Computing relative t0 wrt to chamber average" << endl;
00419   //Compute the reference for each chamber
00420   map<DTChamberId,double> sumT0ByChamber;
00421   map<DTChamberId,int> countT0ByChamber;
00422   for(DTT0::const_iterator tzero = t0s->begin();
00423       tzero != t0s->end(); tzero++) {
00424 // @@@ NEW DTT0 FORMAT
00425 //    DTChamberId chamberId((*tzero).first.wheelId,
00426 //                        (*tzero).first.stationId,
00427 //                        (*tzero).first.sectorId);
00428 //    sumT0ByChamber[chamberId] = sumT0ByChamber[chamberId] + (*tzero).second.t0mean;
00429     int channelId = tzero->channelId;
00430     if ( channelId == 0 ) continue;
00431     DTWireId wireId(channelId);
00432     DTChamberId chamberId(wireId.chamberId());
00433     sumT0ByChamber[chamberId] = sumT0ByChamber[chamberId] + tzero->t0mean;
00434 // @@@ better DTT0 usage
00435     float t0mean_f;
00436     float t0rms_f;
00437     t0s->get(wireId,t0mean_f,t0rms_f,DTTimeUnits::counts);
00438     sumT0ByChamber[chamberId] = sumT0ByChamber[chamberId] + t0mean_f;
00439 // @@@ NEW DTT0 END
00440     countT0ByChamber[chamberId]++;
00441   }
00442 
00443   //Change reference for each wire and store the new t0s in the new map
00444   for(DTT0::const_iterator tzero = t0s->begin();
00445       tzero != t0s->end(); tzero++) {
00446 // @@@ NEW DTT0 FORMAT
00447 //    DTChamberId chamberId((*tzero).first.wheelId,
00448 //                        (*tzero).first.stationId,
00449 //                        (*tzero).first.sectorId);
00450 //    double t0mean = ((*tzero).second.t0mean) - (sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00451 //    double t0rms = (*tzero).second.t0rms;
00452 //    DTWireId wireId((*tzero).first.wheelId,
00453 //                  (*tzero).first.stationId,
00454 //                  (*tzero).first.sectorId,
00455 //                  (*tzero).first.slId,
00456 //                  (*tzero).first.layerId,
00457 //                  (*tzero).first.cellId);
00458     int channelId = tzero->channelId;
00459     if ( channelId == 0 ) continue;
00460     DTWireId wireId(channelId);
00461     DTChamberId chamberId(wireId.chamberId());
00462     double t0mean = (tzero->t0mean) - (sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00463     double t0rms = tzero->t0rms;
00464 // @@@ better DTT0 usage
00465     float t0mean_f;
00466     float t0rms_f;
00467     t0s->get(wireId,t0mean_f,t0rms_f,DTTimeUnits::counts);
00468     t0mean = t0mean_f - (sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00469     t0rms = t0rms_f;
00470 // @@@ NEW DTT0 END
00471     t0sWRTChamber->set(wireId,
00472                        t0mean,
00473                        t0rms,
00474                        DTTimeUnits::counts);
00475     if(debug){
00476       //cout<<"Chamber "<<chamberId<<" has reference "<<(sumT0ByChamber[chamberId]/countT0ByChamber[chamberId]);
00477 //      cout<<"Changing t0 of wire "<<wireId<<" from "<<(*tzero).second.t0mean<<" to "<<t0mean<<endl;
00478       cout<<"Changing t0 of wire "<<wireId<<" from "<<tzero->t0mean<<" to "<<t0mean<<endl;
00479     }
00480   }
00481 
00483   if(debug) 
00484    cout << "[DTT0Calibration]Writing values in DB!" << endl;
00485   // FIXME: to be read from cfg?
00486   string t0Record = "DTT0Rcd";
00487   // Write the t0 map to DB
00488   DTCalibDBUtils::writeToDB(t0Record, t0sWRTChamber);
00489 }
00490 
00491 string DTT0Calibration::getHistoName(const DTWireId& wId) const {
00492   string histoName;
00493   stringstream theStream;
00494   theStream << "Ch_" << wId.wheel() << "_" << wId.station() << "_" << wId.sector()
00495             << "_SL" << wId.superlayer() << "_L" << wId.layer() << "_W"<< wId.wire() <<"_hT0Histo";
00496   theStream >> histoName;
00497   return histoName;
00498 }
00499 
00500 string DTT0Calibration::getHistoName(const DTLayerId& lId) const {
00501   string histoName;
00502   stringstream theStream;
00503   theStream << "Ch_" << lId.wheel() << "_" << lId.station() << "_" << lId.sector()
00504             << "_SL" << lId.superlayer() << "_L" << lId.layer() <<"_hT0Histo";
00505   theStream >> histoName;
00506   return histoName;
00507 }
00508